/home/runner/work/kynema-sgf/kynema-sgf/src/fvm/filter.H Source File

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/fvm/filter.H Source File
Kynema-SGF 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
4#include "src/fvm/fvm_utils.H"
5#include "AMReX_REAL.H"
6
7using namespace amrex::literals;
8
9namespace kynema_sgf::fvm {
10
14template <typename FTypeIn, typename FTypeOut>
15struct Filter
16{
21 Filter(FTypeOut& filterphi, const FTypeIn& phi)
22 : m_filterphi(filterphi), m_phi(phi)
23 {
24 AMREX_ALWAYS_ASSERT(filterphi.num_comp() == phi.num_comp());
25 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
26 }
27
28 template <typename Stencil>
29 void apply(const int lev) const
30 {
31 const int ncomp = m_phi.num_comp();
32 const auto& geom = m_phi.repo().mesh().Geom(lev);
33 const auto dlo = geom.Domain().smallEnd();
34 const auto dhi = geom.Domain().bigEnd();
35
36 auto& out_mf = m_filterphi(lev);
37 auto const& in_arrs = m_phi(lev).const_arrays();
38 auto const& out_arrs = out_mf.arrays();
39
40 amrex::ParallelFor(
41 out_mf, amrex::IntVect(0),
42 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
43 if (!Stencil::applies_to(i, j, k, dlo, dhi)) {
44 return;
45 }
46 auto const& phi_arr = in_arrs[nbx];
47 auto const& filterphi_arr = out_arrs[nbx];
48 for (int icomp = 0; icomp < ncomp; icomp++) {
49 const amrex::Real filx =
50 Stencil::f00 * phi_arr(i + 1, j, k, icomp) +
51 Stencil::f01 * phi_arr(i, j, k, icomp) +
52 Stencil::f02 * phi_arr(i - 1, j, k, icomp);
53 const amrex::Real fily =
54 Stencil::f10 * phi_arr(i, j + 1, k, icomp) +
55 Stencil::f11 * phi_arr(i, j, k, icomp) +
56 Stencil::f12 * phi_arr(i, j - 1, k, icomp);
57 const amrex::Real filz =
58 Stencil::f20 * phi_arr(i, j, k + 1, icomp) +
59 Stencil::f21 * phi_arr(i, j, k, icomp) +
60 Stencil::f22 * phi_arr(i, j, k - 1, icomp);
61 filterphi_arr(i, j, k, icomp) =
62 1.0_rt / 3.0_rt * (filx + fily + filz);
63 }
64 });
65 }
66
67 FTypeOut& m_filterphi;
68 const FTypeIn& m_phi;
69};
70
77template <typename FTypeIn, typename FTypeOut>
78void filter(FTypeOut& filterphi, const FTypeIn& phi)
79{
80 BL_PROFILE("kynema-sgf::fvm::filter");
81 Filter<FTypeIn, FTypeOut> filtering(filterphi, phi);
82 impl::apply(filtering, phi);
83}
84
90template <typename FType>
91std::unique_ptr<ScratchField> filter(const FType& phi)
92{
93 const std::string gname = phi.name() + "_filter";
94 auto filterphi = phi.repo().create_scratch_field(gname, phi.num_comp());
95 filter(*filterphi, phi);
96 return filterphi;
97}
98
99} // namespace kynema_sgf::fvm
100
101#endif /* FILTER_H */
void filter(FTypeOut &filterphi, const FTypeIn &phi)
Definition filter.H:78
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:20
Definition SchemeTraits.H:6
Definition filter.H:16
const FTypeIn & m_phi
Definition filter.H:68
FTypeOut & m_filterphi
Definition filter.H:67
void apply(const int lev) const
Definition filter.H:29
Filter(FTypeOut &filterphi, const FTypeIn &phi)
Definition filter.H:21