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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/nonLinearSum.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
nonLinearSum.H
Go to the documentation of this file.
1#ifndef NONLINEARSUM_H
2#define NONLINEARSUM_H
3
5
6namespace amr_wind::fvm {
7
12template <typename FTypeIn, typename FTypeOut>
14{
15 NonLinearSum(FTypeOut& strphi, const FTypeIn& phi)
16 : m_strphi(strphi), m_phi(phi)
17 {
18 AMREX_ALWAYS_ASSERT(AMREX_SPACEDIM == m_phi.num_comp());
19 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
20 }
21
22 template <typename Stencil>
23 void apply(const int lev, const amrex::MFIter& mfi) const
24 {
25 const auto& geom = m_phi.repo().mesh().Geom(lev);
26 const auto& idx = geom.InvCellSizeArray();
27 const auto& strphi = m_strphi(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 amrex::Real cp1, c, cm1, ux, uy, uz, vx, vy, vz, wx, wy, wz;
39 cp1 = Stencil::c00;
40 c = Stencil::c01;
41 cm1 = Stencil::c02;
42
43 ux = (cp1 * phi(i + 1, j, k, 0) + c * phi(i, j, k, 0) +
44 cm1 * phi(i - 1, j, k, 0)) *
45 idx[0];
46 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
47 cm1 * phi(i - 1, j, k, 1)) *
48 idx[0];
49 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
50 cm1 * phi(i - 1, j, k, 2)) *
51 idx[0];
52
53 cp1 = Stencil::c10;
54 c = Stencil::c11;
55 cm1 = Stencil::c12;
56
57 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
58 cm1 * phi(i, j - 1, k, 0)) *
59 idx[1];
60 vy = (cp1 * phi(i, j + 1, k, 1) + c * phi(i, j, k, 1) +
61 cm1 * phi(i, j - 1, k, 1)) *
62 idx[1];
63 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
64 cm1 * phi(i, j - 1, k, 2)) *
65 idx[1];
66
67 cp1 = Stencil::c20;
68 c = Stencil::c21;
69 cm1 = Stencil::c22;
70
71 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
72 cm1 * phi(i, j, k - 1, 0)) *
73 idx[2];
74 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
75 cm1 * phi(i, j, k - 1, 1)) *
76 idx[2];
77 wz = (cp1 * phi(i, j, k + 1, 2) + c * phi(i, j, k, 2) +
78 cm1 * phi(i, j, k - 1, 2)) *
79 idx[2];
80
81 // N11
82 strphi(i, j, k, 0) = (ux * ux + uy * vx + uz * wx) -
83 (ux * ux + uy * uy + uz * uz) +
84 3 * (ux * ux + vx * vx + wx * wx) +
85 (ux * ux + vx * uy + wx * uz);
86 // N12
87 strphi(i, j, k, 1) = (ux * uy + uy * vy + uz * wy) -
88 (ux * vx + uy * vy + uz * vz) +
89 3 * (ux * uy + vx * vy + wx * wy) +
90 (ux * vx + vx * vy + wx * vz);
91 // N13
92 strphi(i, j, k, 2) = (ux * uz + uy * vz + uz * wz) -
93 (ux * wx + uy * wy + uz * wz) +
94 3 * (ux * uz + vx * vz + wx * wz) +
95 (ux * wx + vx * wy + wx * wz);
96 // N21
97 strphi(i, j, k, 3) = (vx * ux + vy * vx + vz * wx) -
98 (vx * ux + vy * uy + vz * uz) +
99 3 * (uy * ux + vy * vx + wy * wx) +
100 (uy * ux + vy * uy + wy * uz);
101 // N22
102 strphi(i, j, k, 4) = (vx * uy + vy * vy + vz * wy) -
103 (vx * vx + vy * vy + vz * vz) +
104 3 * (uy * vx + vy * vy + wy * wy) +
105 (uy * vx + vy * vy + wy * vz);
106 // N23
107 strphi(i, j, k, 5) = (vx * uz + vy * vz + vz * wz) -
108 (vx * wx + vy * wy + vz * wz) +
109 3 * (uy * wx + vy * vz + wy * wz) +
110 (uy * wx + vy * wy + wy * wz);
111 // N31
112 strphi(i, j, k, 6) = (wx * ux + wy * vx + wz * wx) -
113 (wx * ux + wy * uy + wz * uz) +
114 3 * (uz * ux + vz * vx + wz * wx) +
115 (uz * ux + vz * uy + wz * uz);
116 // N32
117 strphi(i, j, k, 7) = (wx * uy + wy * vy + wz * wy) -
118 (wx * vx + wy * vy + wz * vz) +
119 3 * (uz * uy + vz * vy + wz * wy) +
120 (uz * vx + vz * vy + wz * vz);
121 // N33
122 strphi(i, j, k, 8) = (wx * uz + wy * vz + wz * wz) -
123 (wx * wx + wy * wy + wz * wz) +
124 3 * (uz * uz + vz * vz + wz * wz) +
125 (uz * wx + vz * wy + wz * wz);
126 });
127 }
128
129 FTypeOut& m_strphi;
130 const FTypeIn& m_phi;
131};
132
139template <typename FTypeIn, typename FTypeOut>
140inline void nonlinearsum(FTypeOut& strphi, const FTypeIn& phi)
141{
142 BL_PROFILE("amr-wind::fvm::nonlinearsum");
143 NonLinearSum<FTypeIn, FTypeOut> str(strphi, phi);
144 impl::apply(str, phi);
145}
146
152template <typename FType>
153inline std::unique_ptr<ScratchField> nonlinearsum(const FType& phi)
154{
155 const std::string gname = phi.name() + "_nonlinearsum";
156 auto strphi = phi.repo().create_scratch_field(gname, 9);
157 nonlinearsum(*strphi, phi);
158 return strphi;
159}
160
161} // namespace amr_wind::fvm
162
163#endif /* NONLINEARSUM_H */
void nonlinearsum(FTypeOut &strphi, const FTypeIn &phi)
Definition nonLinearSum.H:140
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition nonLinearSum.H:14
void apply(const int lev, const amrex::MFIter &mfi) const
Definition nonLinearSum.H:23
const FTypeIn & m_phi
Definition nonLinearSum.H:130
FTypeOut & m_strphi
Definition nonLinearSum.H:129
NonLinearSum(FTypeOut &strphi, const FTypeIn &phi)
Definition nonLinearSum.H:15