/home/runner/work/amr-wind/amr-wind/amr-wind/fvm/filter.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/filter.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
filter.H
Go to the documentation of this file.
1#ifndef FILTER_H
2#define FILTER_H
3
5#include "AMReX_Array4.H"
6#include "AMReX_Geometry.H"
7#include "AMReX_MFIter.H"
8#include "AMReX_REAL.H"
9
10using namespace amrex::literals;
11
12namespace amr_wind::fvm {
13
17template <typename FTypeIn, typename FTypeOut>
18struct Filter
19{
24 Filter(FTypeOut& filterphi, const FTypeIn& phi)
25 : m_filterphi(filterphi), m_phi(phi)
26 {
27 AMREX_ALWAYS_ASSERT(filterphi.num_comp() == phi.num_comp());
28 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
29 }
30
31 template <typename Stencil>
32 void apply(const int lev, const amrex::MFIter& mfi) const
33 {
34 const int ncomp = m_phi.num_comp();
35 const auto& geom = m_phi.repo().mesh().Geom(lev);
36 const auto& filterphi_arr = m_filterphi(lev).array(mfi);
37 const auto& phi_arr = m_phi(lev).const_array(mfi);
38
39 const auto& bx_in = mfi.tilebox();
40 const auto& bx = Stencil::box(bx_in, geom);
41 if (bx.isEmpty()) {
42 return;
43 }
44
45 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
46 for (int icomp = 0; icomp < ncomp; icomp++) {
47 amrex::Real fp1 = Stencil::f00;
48 amrex::Real f = Stencil::f01;
49 amrex::Real fm1 = Stencil::f02;
50 const amrex::Real filx =
51 ((fp1 * phi_arr(i + 1, j, k, icomp)) +
52 (f * phi_arr(i, j, k, icomp)) +
53 (fm1 * phi_arr(i - 1, j, k, icomp)));
54
55 fp1 = Stencil::f10;
56 f = Stencil::f11;
57 fm1 = Stencil::f12;
58 const amrex::Real fily =
59 ((fp1 * phi_arr(i, j + 1, k, icomp)) +
60 (f * phi_arr(i, j, k, icomp)) +
61 (fm1 * phi_arr(i, j - 1, k, icomp)));
62
63 fp1 = Stencil::f20;
64 f = Stencil::f21;
65 fm1 = Stencil::f22;
66 const amrex::Real filz =
67 ((fp1 * phi_arr(i, j, k + 1, icomp)) +
68 (f * phi_arr(i, j, k, icomp)) +
69 (fm1 * phi_arr(i, j, k - 1, icomp)));
70 filterphi_arr(i, j, k, icomp) =
71 1.0_rt / 3.0_rt * (filx + fily + filz);
72 }
73 });
74 }
75
76 FTypeOut& m_filterphi;
77 const FTypeIn& m_phi;
78};
79
86template <typename FTypeIn, typename FTypeOut>
87void filter(FTypeOut& filterphi, const FTypeIn& phi)
88{
89 BL_PROFILE("amr-wind::fvm::filter");
90 Filter<FTypeIn, FTypeOut> filtering(filterphi, phi);
91 impl::apply(filtering, phi);
92}
93
99template <typename FType>
100std::unique_ptr<ScratchField> filter(const FType& phi)
101{
102 const std::string gname = phi.name() + "_filter";
103 auto filterphi = phi.repo().create_scratch_field(gname, phi.num_comp());
104 filter(*filterphi, phi);
105 return filterphi;
106}
107
108} // namespace amr_wind::fvm
109
110#endif /* FILTER_H */
void filter(FTypeOut &filterphi, const FTypeIn &phi)
Definition filter.H:87
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition filter.H:19
FTypeOut & m_filterphi
Definition filter.H:76
Filter(FTypeOut &filterphi, const FTypeIn &phi)
Definition filter.H:24
void apply(const int lev, const amrex::MFIter &mfi) const
Definition filter.H:32
const FTypeIn & m_phi
Definition filter.H:77