/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(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
45 for (int icomp = 0; icomp < ncomp; icomp++) {
46 amrex::Real cp1 = Stencil::c00;
47 amrex::Real c = Stencil::c01;
48 amrex::Real cm1 = Stencil::c02;
49 gradphi_arr(i, j, k, (icomp * AMREX_SPACEDIM) + 0) =
50 (cp1 * phi_arr(i + 1, j, k, icomp) +
51 c * phi_arr(i, j, k, icomp) +
52 cm1 * phi_arr(i - 1, j, k, icomp)) *
53 idx[0];
54
55 cp1 = Stencil::c10;
56 c = Stencil::c11;
57 cm1 = Stencil::c12;
58 gradphi_arr(i, j, k, (icomp * AMREX_SPACEDIM) + 1) =
59 (cp1 * phi_arr(i, j + 1, k, icomp) +
60 c * phi_arr(i, j, k, icomp) +
61 cm1 * phi_arr(i, j - 1, k, icomp)) *
62 idx[1];
63
64 cp1 = Stencil::c20;
65 c = Stencil::c21;
66 cm1 = Stencil::c22;
67 gradphi_arr(i, j, k, (icomp * AMREX_SPACEDIM) + 2) =
68 (cp1 * phi_arr(i, j, k + 1, icomp) +
69 c * phi_arr(i, j, k, icomp) +
70 cm1 * phi_arr(i, j, k - 1, icomp)) *
71 idx[2];
72 }
73 });
74 }
75
76 FTypeOut& m_gradphi;
77 const FTypeIn& m_phi;
78};
79
86template <typename FTypeIn, typename FTypeOut>
87void gradient(FTypeOut& gradphi, const FTypeIn& phi)
88{
89 BL_PROFILE("amr-wind::fvm::gradient");
90 Gradient<FTypeIn, FTypeOut> grad(gradphi, phi);
91 impl::apply(grad, phi);
92}
93
99template <typename FType>
100std::unique_ptr<ScratchField> gradient(const FType& phi)
101{
102 const std::string gname = phi.name() + "_gradient";
103 auto gradphi =
104 phi.repo().create_scratch_field(gname, phi.num_comp() * AMREX_SPACEDIM);
105 gradient(*gradphi, phi);
106 return gradphi;
107}
108
109} // namespace amr_wind::fvm
110
111#endif /* GRADIENT_H */
void gradient(FTypeOut &gradphi, const FTypeIn &phi)
Definition gradient.H:87
AMREX_INLINE 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:77
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:76