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

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/fvm/divergence.H Source File
Kynema-SGF API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
divergence.H
Go to the documentation of this file.
1#ifndef DIVERGENCE_H
2#define DIVERGENCE_H
3
4#include "src/fvm/fvm_utils.H"
5
6namespace kynema_sgf::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
14 Divergence(FTypeOut& divphi, const FTypeIn& phi)
15 : m_divphi(divphi), m_phi(phi)
16 {
17 AMREX_ALWAYS_ASSERT(
18 AMREX_SPACEDIM * divphi.num_comp() == phi.num_comp());
19 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
20 }
21
22 template <typename Stencil>
23 void apply(const int lev) const
24 {
25 const int ncomp = m_divphi.num_comp();
26 const auto& geom = m_phi.repo().mesh().Geom(lev);
27 const auto& idx = geom.InvCellSizeArray();
28 const auto dlo = geom.Domain().smallEnd();
29 const auto dhi = geom.Domain().bigEnd();
30
31 auto& out_mf = m_divphi(lev);
32 auto const& in_arrs = m_phi(lev).const_arrays();
33 auto const& out_arrs = out_mf.arrays();
34
35 amrex::ParallelFor(
36 out_mf, amrex::IntVect(0),
37 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
38 if (!Stencil::applies_to(i, j, k, dlo, dhi)) {
39 return;
40 }
41 auto const& phi_arr = in_arrs[nbx];
42 auto const& divphi_arr = out_arrs[nbx];
43 for (int icomp = 0; icomp < ncomp; icomp++) {
44 divphi_arr(i, j, k, icomp) =
45 (Stencil::c00 *
46 phi_arr(
47 i + 1, j, k, (icomp * AMREX_SPACEDIM) + 0) +
48 Stencil::c01 *
49 phi_arr(i, j, k, (icomp * AMREX_SPACEDIM) + 0) +
50 Stencil::c02 *
51 phi_arr(
52 i - 1, j, k, (icomp * AMREX_SPACEDIM) + 0)) *
53 idx[0];
54
55 divphi_arr(i, j, k, icomp) +=
56 (Stencil::c10 *
57 phi_arr(
58 i, j + 1, k, (icomp * AMREX_SPACEDIM) + 1) +
59 Stencil::c11 *
60 phi_arr(i, j, k, (icomp * AMREX_SPACEDIM) + 1) +
61 Stencil::c12 *
62 phi_arr(
63 i, j - 1, k, (icomp * AMREX_SPACEDIM) + 1)) *
64 idx[1];
65
66 divphi_arr(i, j, k, icomp) +=
67 (Stencil::c20 *
68 phi_arr(
69 i, j, k + 1, (icomp * AMREX_SPACEDIM) + 2) +
70 Stencil::c21 *
71 phi_arr(i, j, k, (icomp * AMREX_SPACEDIM) + 2) +
72 Stencil::c22 *
73 phi_arr(
74 i, j, k - 1, (icomp * AMREX_SPACEDIM) + 2)) *
75 idx[2];
76 }
77 });
78 }
79
80 FTypeOut& m_divphi;
81 const FTypeIn& m_phi;
82};
83
87template <typename FTypeIn, typename FTypeOut>
88void divergence(FTypeOut& divphi, const FTypeIn& phi)
89{
90 BL_PROFILE("kynema-sgf::fvm::divergence");
91 Divergence<FTypeIn, FTypeOut> grad(divphi, phi);
92 impl::apply(grad, phi);
93}
94
98template <typename FType>
99std::unique_ptr<ScratchField> divergence(const FType& phi)
100{
101 AMREX_ALWAYS_ASSERT(phi.num_comp() >= AMREX_SPACEDIM);
102 AMREX_ALWAYS_ASSERT(phi.num_comp() % AMREX_SPACEDIM == 0);
103 const std::string gname = phi.name() + "_divergence";
104 auto divphi =
105 phi.repo().create_scratch_field(gname, phi.num_comp() / AMREX_SPACEDIM);
106 divergence(*divphi, phi);
107 return divphi;
108}
109
110} // namespace kynema_sgf::fvm
111
112#endif /* DIVERGENCE_H */
void divergence(FTypeOut &divphi, const FTypeIn &phi)
Definition divergence.H:88
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:20
Definition SchemeTraits.H:6
Definition divergence.H:13
FTypeOut & m_divphi
Definition divergence.H:80
const FTypeIn & m_phi
Definition divergence.H:81
void apply(const int lev) const
Definition divergence.H:23
Divergence(FTypeOut &divphi, const FTypeIn &phi)
Definition divergence.H:14