/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
6namespace amr_wind::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
14 Curvature(FTypeOut& curphi, const FTypeIn& phi)
15 : m_curphi(curphi), m_phi(phi)
16 {
17 AMREX_ALWAYS_ASSERT(m_phi.num_comp() == m_curphi.num_comp());
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 auto& geom = m_phi.repo().mesh().Geom(lev);
25 const auto& idx = geom.InvCellSizeArray();
26 const auto& curphi = m_curphi(lev).array(mfi);
27 const auto& phi = m_phi(lev).const_array(mfi);
28
29 const auto& bx_in = mfi.tilebox();
30 const auto& bx = Stencil::box(bx_in, geom);
31 if (bx.isEmpty()) {
32 return;
33 }
34
35 amrex::ParallelFor(
36 bx, m_phi.num_comp(),
37 [=] AMREX_GPU_DEVICE(int i, int j, int k, int n) noexcept {
38 amrex::Real cp1, c, cm1;
39 amrex::Real sp1, s, sm1;
40
41 sp1 = Stencil::s00;
42 s = Stencil::s01;
43 sm1 = Stencil::s02;
44 const amrex::Real phixx =
45 (sp1 * phi(i + 1, j, k, n) + s * phi(i, j, k, n) +
46 sm1 * phi(i - 1, j, k, n)) *
47 idx[0] * idx[0];
48
49 sp1 = Stencil::s10;
50 s = Stencil::s11;
51 sm1 = Stencil::s12;
52 const amrex::Real phiyy =
53 (sp1 * phi(i, j + 1, k, n) + s * phi(i, j, k, n) +
54 sm1 * phi(i, j - 1, k, n)) *
55 idx[1] * idx[1];
56
57 sp1 = Stencil::s20;
58 s = Stencil::s21;
59 sm1 = Stencil::s22;
60 const amrex::Real phizz =
61 (sp1 * phi(i, j, k + 1, n) + s * phi(i, j, k, n) +
62 sm1 * phi(i, j, k - 1, n)) *
63 idx[2] * idx[2];
64
65 cp1 = Stencil::c20;
66 c = Stencil::c21;
67 cm1 = Stencil::c22;
68 const amrex::Real phiz =
69 (cp1 * phi(i, j, k + 1, n) + c * phi(i, j, k, n) +
70 cm1 * phi(i, j, k - 1, n)) *
71 idx[2];
72 const amrex::Real phiz_ip1 =
73 (cp1 * phi(i + 1, j, k + 1, n) + c * phi(i + 1, j, k, n) +
74 cm1 * phi(i + 1, j, k - 1, n)) *
75 idx[2];
76 const amrex::Real phiz_im1 =
77 (cp1 * phi(i - 1, j, k + 1, n) + c * phi(i - 1, j, k, n) +
78 cm1 * phi(i - 1, j, k - 1, n)) *
79 idx[2];
80 const amrex::Real phiz_jp1 =
81 (cp1 * phi(i, j + 1, k + 1, n) + c * phi(i, j + 1, k, n) +
82 cm1 * phi(i, j + 1, k - 1, n)) *
83 idx[2];
84 const amrex::Real phiz_jm1 =
85 (cp1 * phi(i, j - 1, k + 1, n) + c * phi(i, j - 1, k, n) +
86 cm1 * phi(i, j - 1, k - 1, n)) *
87 idx[2];
88
89 cp1 = Stencil::c10;
90 c = Stencil::c11;
91 cm1 = Stencil::c12;
92 const amrex::Real phiy =
93 (cp1 * phi(i, j + 1, k, n) + c * phi(i, j, k, n) +
94 cm1 * phi(i, j - 1, k, n)) *
95 idx[1];
96 const amrex::Real phiy_ip1 =
97 (cp1 * phi(i + 1, j + 1, k, n) + c * phi(i + 1, j, k, n) +
98 cm1 * phi(i + 1, j - 1, k, n)) *
99 idx[1];
100 const amrex::Real phiy_im1 =
101 (cp1 * phi(i - 1, j + 1, k, n) + c * phi(i - 1, j, k, n) +
102 cm1 * phi(i - 1, j - 1, k, n)) *
103 idx[1];
104 const amrex::Real phiyz =
105 (cp1 * phiz_jp1 + c * phiz + cm1 * phiz_jm1) * idx[1];
106
107 cp1 = Stencil::c00;
108 c = Stencil::c01;
109 cm1 = Stencil::c02;
110 const amrex::Real phix =
111 (cp1 * phi(i + 1, j, k, n) + c * phi(i, j, k, n) +
112 cm1 * phi(i - 1, j, k, n)) *
113 idx[0];
114 const amrex::Real phixy =
115 (cp1 * phiy_ip1 + c * phiy + cm1 * phiy_im1) * idx[0];
116 const amrex::Real phixz =
117 (cp1 * phiz_ip1 + c * phiz + cm1 * phiz_im1) * idx[0];
118
119 curphi(i, j, k, n) =
120 -(phix * phix * phiyy - 2. * phix * phiy * phixy +
121 phiy * phiy * phixx + phix * phix * phizz -
122 2. * phix * phiz * phixz + phiz * phiz * phixx +
123 phiy * phiy * phizz - 2. * phiy * phiz * phiyz +
124 phiz * phiz * phiyy) /
125 std::pow(phix * phix + phiy * phiy + phiz * phiz, 1.5);
126 });
127 }
128
129 FTypeOut& m_curphi;
130 const FTypeIn& m_phi;
131};
132
136template <typename FTypeIn, typename FTypeOut>
137inline void curvature(FTypeOut& curphi, const FTypeIn& phi)
138{
139 BL_PROFILE("amr-wind::fvm::curvature");
140 Curvature<FTypeIn, FTypeOut> cur(curphi, phi);
141 impl::apply(cur, phi);
142}
143
147template <typename FType>
148inline std::unique_ptr<ScratchField> curvature(const FType& phi)
149{
150 const std::string gname = phi.name() + "_curvature";
151 auto curphi = phi.repo().create_scratch_field(gname, phi.num_comp());
152 curvature(*curphi, phi);
153 return curphi;
154}
155
156} // namespace amr_wind::fvm
157
158#endif /* CURVATURE_H */
void curvature(FTypeOut &curphi, const FTypeIn &phi)
Definition curvature.H:137
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition curvature.H:13
FTypeOut & m_curphi
Definition curvature.H:129
Curvature(FTypeOut &curphi, const FTypeIn &phi)
Definition curvature.H:14
void apply(const int lev, const amrex::MFIter &mfi) const
Definition curvature.H:22
const FTypeIn & m_phi
Definition curvature.H:130