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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/strainrate.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
strainrate.H
Go to the documentation of this file.
1#ifndef STRAINRATE_H
2#define STRAINRATE_H
3
5
6namespace amr_wind::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
14 StrainRate(FTypeOut& strphi, const FTypeIn& phi)
15 : m_strphi(strphi), m_phi(phi)
16 {
17 AMREX_ALWAYS_ASSERT(AMREX_SPACEDIM == m_phi.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& strphi = m_strphi(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, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
37 amrex::Real cp1, c, cm1, ux, uy, uz, vx, vy, vz, wx, wy, wz;
38 cp1 = Stencil::c00;
39 c = Stencil::c01;
40 cm1 = Stencil::c02;
41
42 ux = (cp1 * phi(i + 1, j, k, 0) + c * phi(i, j, k, 0) +
43 cm1 * phi(i - 1, j, k, 0)) *
44 idx[0];
45 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
46 cm1 * phi(i - 1, j, k, 1)) *
47 idx[0];
48 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
49 cm1 * phi(i - 1, j, k, 2)) *
50 idx[0];
51
52 cp1 = Stencil::c10;
53 c = Stencil::c11;
54 cm1 = Stencil::c12;
55
56 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
57 cm1 * phi(i, j - 1, k, 0)) *
58 idx[1];
59 vy = (cp1 * phi(i, j + 1, k, 1) + c * phi(i, j, k, 1) +
60 cm1 * phi(i, j - 1, k, 1)) *
61 idx[1];
62 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
63 cm1 * phi(i, j - 1, k, 2)) *
64 idx[1];
65
66 cp1 = Stencil::c20;
67 c = Stencil::c21;
68 cm1 = Stencil::c22;
69
70 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
71 cm1 * phi(i, j, k - 1, 0)) *
72 idx[2];
73 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
74 cm1 * phi(i, j, k - 1, 1)) *
75 idx[2];
76 wz = (cp1 * phi(i, j, k + 1, 2) + c * phi(i, j, k, 2) +
77 cm1 * phi(i, j, k - 1, 2)) *
78 idx[2];
79
80 strphi(i, j, k) = sqrt(
81 2.0 * std::pow(ux, 2) + 2.0 * std::pow(vy, 2) +
82 2.0 * std::pow(wz, 2) + std::pow(uy + vx, 2) +
83 std::pow(vz + wy, 2) + std::pow(wx + uz, 2));
84 });
85 }
86
87 FTypeOut& m_strphi;
88 const FTypeIn& m_phi;
89};
90
97template <typename FTypeIn, typename FTypeOut>
98inline void strainrate(FTypeOut& strphi, const FTypeIn& phi)
99{
100 BL_PROFILE("amr-wind::fvm::strainrate");
101 StrainRate<FTypeIn, FTypeOut> str(strphi, phi);
102 impl::apply(str, phi);
103}
104
110template <typename FType>
111inline std::unique_ptr<ScratchField> strainrate(const FType& phi)
112{
113 const std::string gname = phi.name() + "_strainrate";
114 auto strphi = phi.repo().create_scratch_field(gname, 1);
115 strainrate(*strphi, phi);
116 return strphi;
117}
118
119} // namespace amr_wind::fvm
120
121#endif /* STRAINRATE_H */
void strainrate(FTypeOut &strphi, const FTypeIn &phi)
Definition strainrate.H:98
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition strainrate.H:13
const FTypeIn & m_phi
Definition strainrate.H:88
void apply(const int lev, const amrex::MFIter &mfi) const
Definition strainrate.H:22
FTypeOut & m_strphi
Definition strainrate.H:87
StrainRate(FTypeOut &strphi, const FTypeIn &phi)
Definition strainrate.H:14