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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/divergence.H Source File
AMR-Wind 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
5#include "AMReX_Array4.H"
6#include "AMReX_Geometry.H"
7#include "AMReX_MFIter.H"
8
9namespace amr_wind::fvm {
10
14template <typename FTypeIn, typename FTypeOut>
16{
17 Divergence(FTypeOut& divphi, const FTypeIn& phi)
18 : m_divphi(divphi), m_phi(phi)
19 {
20 AMREX_ALWAYS_ASSERT(
21 AMREX_SPACEDIM * divphi.num_comp() == phi.num_comp());
22 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
23 }
24
25 template <typename Stencil>
26 void apply(const int lev, const amrex::MFIter& mfi) const
27 {
28 const int ncomp = m_divphi.num_comp();
29 const auto& geom = m_phi.repo().mesh().Geom(lev);
30 const auto& idx = geom.InvCellSizeArray();
31 const auto& divphi_arr = m_divphi(lev).array(mfi);
32 const auto& phi_arr = m_phi(lev).const_array(mfi);
33
34 const auto& bx_in = mfi.tilebox();
35 const auto& bx = Stencil::box(bx_in, geom);
36 if (bx.isEmpty()) {
37 return;
38 }
39
40 amrex::ParallelFor(
41 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
42 for (int icomp = 0; icomp < ncomp; icomp++) {
43 amrex::Real cp1 = Stencil::c00;
44 amrex::Real c = Stencil::c01;
45 amrex::Real cm1 = Stencil::c02;
46 divphi_arr(i, j, k, icomp) =
47 (cp1 *
48 phi_arr(i + 1, j, k, icomp * AMREX_SPACEDIM + 0) +
49 c * phi_arr(i, j, k, icomp * AMREX_SPACEDIM + 0) +
50 cm1 *
51 phi_arr(i - 1, j, k, icomp * AMREX_SPACEDIM + 0)) *
52 idx[0];
53
54 cp1 = Stencil::c10;
55 c = Stencil::c11;
56 cm1 = Stencil::c12;
57 divphi_arr(i, j, k, icomp) +=
58 (cp1 *
59 phi_arr(i, j + 1, k, icomp * AMREX_SPACEDIM + 1) +
60 c * phi_arr(i, j, k, icomp * AMREX_SPACEDIM + 1) +
61 cm1 *
62 phi_arr(i, j - 1, k, icomp * AMREX_SPACEDIM + 1)) *
63 idx[1];
64
65 cp1 = Stencil::c20;
66 c = Stencil::c21;
67 cm1 = Stencil::c22;
68 divphi_arr(i, j, k, icomp) +=
69 (cp1 *
70 phi_arr(i, j, k + 1, icomp * AMREX_SPACEDIM + 2) +
71 c * phi_arr(i, j, k, icomp * AMREX_SPACEDIM + 2) +
72 cm1 *
73 phi_arr(i, j, k - 1, icomp * AMREX_SPACEDIM + 2)) *
74 idx[2];
75 }
76 });
77 }
78
79 FTypeOut& m_divphi;
80 const FTypeIn& m_phi;
81};
82
86template <typename FTypeIn, typename FTypeOut>
87inline void divergence(FTypeOut& divphi, const FTypeIn& phi)
88{
89 BL_PROFILE("amr-wind::fvm::divergence");
90 Divergence<FTypeIn, FTypeOut> grad(divphi, phi);
91 impl::apply(grad, phi);
92}
93
97template <typename FType>
98inline std::unique_ptr<ScratchField> divergence(const FType& phi)
99{
100 AMREX_ALWAYS_ASSERT(phi.num_comp() >= AMREX_SPACEDIM);
101 AMREX_ALWAYS_ASSERT(phi.num_comp() % AMREX_SPACEDIM == 0);
102 const std::string gname = phi.name() + "_divergence";
103 auto divphi =
104 phi.repo().create_scratch_field(gname, phi.num_comp() / AMREX_SPACEDIM);
105 divergence(*divphi, phi);
106 return divphi;
107}
108
109} // namespace amr_wind::fvm
110
111#endif /* DIVERGENCE_H */
void divergence(FTypeOut &divphi, const FTypeIn &phi)
Definition divergence.H:87
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition divergence.H:16
Divergence(FTypeOut &divphi, const FTypeIn &phi)
Definition divergence.H:17
FTypeOut & m_divphi
Definition divergence.H:79
void apply(const int lev, const amrex::MFIter &mfi) const
Definition divergence.H:26
const FTypeIn & m_phi
Definition divergence.H:80