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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/curvature.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
curvature.H
Go to the documentation of this file.
1#ifndef CURVATURE_H
2#define CURVATURE_H
3
5#include "AMReX_REAL.H"
6
7using namespace amrex::literals;
8
9namespace amr_wind::fvm {
10
14template <typename FTypeIn, typename FTypeOut>
16{
17 Curvature(FTypeOut& curphi, const FTypeIn& phi)
18 : m_curphi(curphi), m_phi(phi)
19 {
20 AMREX_ALWAYS_ASSERT(m_phi.num_comp() == m_curphi.num_comp());
21 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
22 }
23
24 template <typename Stencil>
25 void apply(const int lev, const amrex::MFIter& mfi) const
26 {
27 const auto& geom = m_phi.repo().mesh().Geom(lev);
28 const auto& idx = geom.InvCellSizeArray();
29 const auto& curphi = m_curphi(lev).array(mfi);
30 const auto& phi = m_phi(lev).const_array(mfi);
31
32 const auto& bx_in = mfi.tilebox();
33 const auto& bx = Stencil::box(bx_in, geom);
34 if (bx.isEmpty()) {
35 return;
36 }
37
38 amrex::ParallelFor(
39 bx, m_phi.num_comp(),
40 [=] AMREX_GPU_DEVICE(int i, int j, int k, int n) noexcept {
41 amrex::Real cp1, c, cm1;
42 amrex::Real sp1, s, sm1;
43
44 sp1 = Stencil::s00;
45 s = Stencil::s01;
46 sm1 = Stencil::s02;
47 const amrex::Real phixx =
48 (sp1 * phi(i + 1, j, k, n) + s * phi(i, j, k, n) +
49 sm1 * phi(i - 1, j, k, n)) *
50 idx[0] * idx[0];
51
52 sp1 = Stencil::s10;
53 s = Stencil::s11;
54 sm1 = Stencil::s12;
55 const amrex::Real phiyy =
56 (sp1 * phi(i, j + 1, k, n) + s * phi(i, j, k, n) +
57 sm1 * phi(i, j - 1, k, n)) *
58 idx[1] * idx[1];
59
60 sp1 = Stencil::s20;
61 s = Stencil::s21;
62 sm1 = Stencil::s22;
63 const amrex::Real phizz =
64 (sp1 * phi(i, j, k + 1, n) + s * phi(i, j, k, n) +
65 sm1 * phi(i, j, k - 1, n)) *
66 idx[2] * idx[2];
67
68 cp1 = Stencil::c20;
69 c = Stencil::c21;
70 cm1 = Stencil::c22;
71 const amrex::Real phiz =
72 (cp1 * phi(i, j, k + 1, n) + c * phi(i, j, k, n) +
73 cm1 * phi(i, j, k - 1, n)) *
74 idx[2];
75 const amrex::Real phiz_ip1 =
76 (cp1 * phi(i + 1, j, k + 1, n) + c * phi(i + 1, j, k, n) +
77 cm1 * phi(i + 1, j, k - 1, n)) *
78 idx[2];
79 const amrex::Real phiz_im1 =
80 (cp1 * phi(i - 1, j, k + 1, n) + c * phi(i - 1, j, k, n) +
81 cm1 * phi(i - 1, j, k - 1, n)) *
82 idx[2];
83 const amrex::Real phiz_jp1 =
84 (cp1 * phi(i, j + 1, k + 1, n) + c * phi(i, j + 1, k, n) +
85 cm1 * phi(i, j + 1, k - 1, n)) *
86 idx[2];
87 const amrex::Real phiz_jm1 =
88 (cp1 * phi(i, j - 1, k + 1, n) + c * phi(i, j - 1, k, n) +
89 cm1 * phi(i, j - 1, k - 1, n)) *
90 idx[2];
91
92 cp1 = Stencil::c10;
93 c = Stencil::c11;
94 cm1 = Stencil::c12;
95 const amrex::Real phiy =
96 (cp1 * phi(i, j + 1, k, n) + c * phi(i, j, k, n) +
97 cm1 * phi(i, j - 1, k, n)) *
98 idx[1];
99 const amrex::Real phiy_ip1 =
100 (cp1 * phi(i + 1, j + 1, k, n) + c * phi(i + 1, j, k, n) +
101 cm1 * phi(i + 1, j - 1, k, n)) *
102 idx[1];
103 const amrex::Real phiy_im1 =
104 (cp1 * phi(i - 1, j + 1, k, n) + c * phi(i - 1, j, k, n) +
105 cm1 * phi(i - 1, j - 1, k, n)) *
106 idx[1];
107 const amrex::Real phiyz =
108 (cp1 * phiz_jp1 + c * phiz + cm1 * phiz_jm1) * idx[1];
109
110 cp1 = Stencil::c00;
111 c = Stencil::c01;
112 cm1 = Stencil::c02;
113 const amrex::Real phix =
114 (cp1 * phi(i + 1, j, k, n) + c * phi(i, j, k, n) +
115 cm1 * phi(i - 1, j, k, n)) *
116 idx[0];
117 const amrex::Real phixy =
118 (cp1 * phiy_ip1 + c * phiy + cm1 * phiy_im1) * idx[0];
119 const amrex::Real phixz =
120 (cp1 * phiz_ip1 + c * phiz + cm1 * phiz_im1) * idx[0];
121
122 curphi(i, j, k, n) =
123 -(phix * phix * phiyy - 2.0_rt * phix * phiy * phixy +
124 phiy * phiy * phixx + phix * phix * phizz -
125 2.0_rt * phix * phiz * phixz + phiz * phiz * phixx +
126 phiy * phiy * phizz - 2.0_rt * phiy * phiz * phiyz +
127 phiz * phiz * phiyy) /
128 std::pow(phix * phix + phiy * phiy + phiz * phiz, 1.5_rt);
129 });
130 }
131
132 FTypeOut& m_curphi;
133 const FTypeIn& m_phi;
134};
135
139template <typename FTypeIn, typename FTypeOut>
140inline void curvature(FTypeOut& curphi, const FTypeIn& phi)
141{
142 BL_PROFILE("amr-wind::fvm::curvature");
143 Curvature<FTypeIn, FTypeOut> cur(curphi, phi);
144 impl::apply(cur, phi);
145}
146
150template <typename FType>
151inline std::unique_ptr<ScratchField> curvature(const FType& phi)
152{
153 const std::string gname = phi.name() + "_curvature";
154 auto curphi = phi.repo().create_scratch_field(gname, phi.num_comp());
155 curvature(*curphi, phi);
156 return curphi;
157}
158
159} // namespace amr_wind::fvm
160
161#endif /* CURVATURE_H */
void curvature(FTypeOut &curphi, const FTypeIn &phi)
Definition curvature.H:140
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition curvature.H:16
FTypeOut & m_curphi
Definition curvature.H:132
Curvature(FTypeOut &curphi, const FTypeIn &phi)
Definition curvature.H:17
void apply(const int lev, const amrex::MFIter &mfi) const
Definition curvature.H:25
const FTypeIn & m_phi
Definition curvature.H:133