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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/laplacian.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
laplacian.H
Go to the documentation of this file.
1#ifndef LAPLACIAN_H
2#define LAPLACIAN_H
3
5
6namespace amr_wind::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
14 Laplacian(FTypeOut& lphi, const FTypeIn& phi) : m_lapphi(lphi), m_phi(phi)
15 {
16 AMREX_ALWAYS_ASSERT(m_lapphi.num_comp() == 1);
17 AMREX_ALWAYS_ASSERT(m_phi.num_comp() == AMREX_SPACEDIM);
18 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
19 }
20
21 template <typename Stencil>
22 void apply(const int lev, const amrex::MFIter& mfi) const
23 {
24 const int ncomp = m_phi.num_comp();
25 const auto& geom = m_phi.repo().mesh().Geom(lev);
26 const auto& idx = geom.InvCellSizeArray();
27 const auto& lapphi = m_lapphi(lev).array(mfi);
28 const auto& phi = m_phi(lev).const_array(mfi);
29
30 const auto& bx_in = mfi.tilebox();
31 const auto& bx = Stencil::box(bx_in, geom);
32 if (bx.isEmpty()) {
33 return;
34 }
35
36 amrex::ParallelFor(
37 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
38 for (int icomp = 0; icomp < ncomp; icomp++) {
39 amrex::Real sp1 = Stencil::s00;
40 amrex::Real s = Stencil::s01;
41 amrex::Real sm1 = Stencil::s02;
42 amrex::Real d2phidx2 =
43 (sp1 * phi(i + 1, j, k, 0) + s * phi(i, j, k, 0) +
44 sm1 * phi(i - 1, j, k, 0)) *
45 idx[0] * idx[0];
46 sp1 = Stencil::s10;
47 s = Stencil::s11;
48 sm1 = Stencil::s12;
49 amrex::Real d2phidy2 =
50 (sp1 * phi(i, j + 1, k, 1) + s * phi(i, j, k, 1) +
51 sm1 * phi(i, j - 1, k, 1)) *
52 idx[1] * idx[1];
53 sp1 = Stencil::s20;
54 s = Stencil::s21;
55 sm1 = Stencil::s22;
56 amrex::Real d2phidz2 =
57 (sp1 * phi(i, j, k + 1, 2) + s * phi(i, j, k, 2) +
58 sm1 * phi(i, j, k - 1, 2)) *
59 idx[2] * idx[2];
60 lapphi(i, j, k) = d2phidx2 + d2phidy2 + d2phidz2;
61 }
62 });
63 }
64
65 FTypeOut& m_lapphi;
66 const FTypeIn& m_phi;
67};
68
75template <typename FTypeIn, typename FTypeOut>
76inline void laplacian(FTypeOut& lapphi, const FTypeIn& phi)
77{
78 BL_PROFILE("amr-wind::fvm::laplacian");
79 Laplacian<FTypeIn, FTypeOut> lap(lapphi, phi);
80 impl::apply(lap, phi);
81}
82
88template <typename FType>
89inline std::unique_ptr<ScratchField> laplacian(const FType& phi)
90{
91 const std::string gname = phi.name() + "_laplacian";
92 auto lapphi = phi.repo().create_scratch_field(gname, 1);
93 laplacian(*lapphi, phi);
94 return lapphi;
95}
96
97} // namespace amr_wind::fvm
98
99#endif /* LAPLACIAN_H */
void laplacian(FTypeOut &lapphi, const FTypeIn &phi)
Definition laplacian.H:76
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition laplacian.H:13
void apply(const int lev, const amrex::MFIter &mfi) const
Definition laplacian.H:22
const FTypeIn & m_phi
Definition laplacian.H:66
FTypeOut & m_lapphi
Definition laplacian.H:65
Laplacian(FTypeOut &lphi, const FTypeIn &phi)
Definition laplacian.H:14