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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/gradient.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
gradient.H
Go to the documentation of this file.
1#ifndef GRADIENT_H
2#define GRADIENT_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{
21 Gradient(FTypeOut& gradphi, const FTypeIn& phi)
22 : m_gradphi(gradphi), m_phi(phi)
23 {
24 AMREX_ALWAYS_ASSERT(
25 gradphi.num_comp() == AMREX_SPACEDIM * phi.num_comp());
26 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
27 }
28
29 template <typename Stencil>
30 void apply(const int lev, const amrex::MFIter& mfi) const
31 {
32 const int ncomp = m_phi.num_comp();
33 const auto& geom = m_phi.repo().mesh().Geom(lev);
34 const auto& idx = geom.InvCellSizeArray();
35 const auto& gradphi_arr = m_gradphi(lev).array(mfi);
36 const auto& phi_arr = m_phi(lev).const_array(mfi);
37
38 const auto& bx_in = mfi.tilebox();
39 const auto& bx = Stencil::box(bx_in, geom);
40 if (bx.isEmpty()) {
41 return;
42 }
43
44 amrex::ParallelFor(
45 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
46 for (int icomp = 0; icomp < ncomp; icomp++) {
47 amrex::Real cp1 = Stencil::c00;
48 amrex::Real c = Stencil::c01;
49 amrex::Real cm1 = Stencil::c02;
50 gradphi_arr(i, j, k, icomp * AMREX_SPACEDIM + 0) =
51 (cp1 * phi_arr(i + 1, j, k, icomp) +
52 c * phi_arr(i, j, k, icomp) +
53 cm1 * phi_arr(i - 1, j, k, icomp)) *
54 idx[0];
55
56 cp1 = Stencil::c10;
57 c = Stencil::c11;
58 cm1 = Stencil::c12;
59 gradphi_arr(i, j, k, icomp * AMREX_SPACEDIM + 1) =
60 (cp1 * phi_arr(i, j + 1, k, icomp) +
61 c * phi_arr(i, j, k, icomp) +
62 cm1 * phi_arr(i, j - 1, k, icomp)) *
63 idx[1];
64
65 cp1 = Stencil::c20;
66 c = Stencil::c21;
67 cm1 = Stencil::c22;
68 gradphi_arr(i, j, k, icomp * AMREX_SPACEDIM + 2) =
69 (cp1 * phi_arr(i, j, k + 1, icomp) +
70 c * phi_arr(i, j, k, icomp) +
71 cm1 * phi_arr(i, j, k - 1, icomp)) *
72 idx[2];
73 }
74 });
75 }
76
77 FTypeOut& m_gradphi;
78 const FTypeIn& m_phi;
79};
80
87template <typename FTypeIn, typename FTypeOut>
88inline void gradient(FTypeOut& gradphi, const FTypeIn& phi)
89{
90 BL_PROFILE("amr-wind::fvm::gradient");
91 Gradient<FTypeIn, FTypeOut> grad(gradphi, phi);
92 impl::apply(grad, phi);
93}
94
100template <typename FType>
101inline std::unique_ptr<ScratchField> gradient(const FType& phi)
102{
103 const std::string gname = phi.name() + "_gradient";
104 auto gradphi =
105 phi.repo().create_scratch_field(gname, phi.num_comp() * AMREX_SPACEDIM);
106 gradient(*gradphi, phi);
107 return gradphi;
108}
109
110} // namespace amr_wind::fvm
111
112#endif /* GRADIENT_H */
void gradient(FTypeOut &gradphi, const FTypeIn &phi)
Definition gradient.H:88
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition gradient.H:16
const FTypeIn & m_phi
Definition gradient.H:78
Gradient(FTypeOut &gradphi, const FTypeIn &phi)
Definition gradient.H:21
void apply(const int lev, const amrex::MFIter &mfi) const
Definition gradient.H:30
FTypeOut & m_gradphi
Definition gradient.H:77