/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#include "AMReX_REAL.H"
6
7using namespace amrex::literals;
8
9namespace amr_wind::fvm {
10
14template <typename FTypeIn, typename FTypeOut>
16{
17 StrainRate(FTypeOut& strphi, const FTypeIn& phi)
18 : m_strphi(strphi), m_phi(phi)
19 {
20 AMREX_ALWAYS_ASSERT(AMREX_SPACEDIM == m_phi.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& strphi = m_strphi(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, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
40 amrex::Real cp1, c, cm1, ux, uy, uz, vx, vy, vz, wx, wy, wz;
41 cp1 = Stencil::c00;
42 c = Stencil::c01;
43 cm1 = Stencil::c02;
44
45 ux = (cp1 * phi(i + 1, j, k, 0) + c * phi(i, j, k, 0) +
46 cm1 * phi(i - 1, j, k, 0)) *
47 idx[0];
48 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
49 cm1 * phi(i - 1, j, k, 1)) *
50 idx[0];
51 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
52 cm1 * phi(i - 1, j, k, 2)) *
53 idx[0];
54
55 cp1 = Stencil::c10;
56 c = Stencil::c11;
57 cm1 = Stencil::c12;
58
59 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
60 cm1 * phi(i, j - 1, k, 0)) *
61 idx[1];
62 vy = (cp1 * phi(i, j + 1, k, 1) + c * phi(i, j, k, 1) +
63 cm1 * phi(i, j - 1, k, 1)) *
64 idx[1];
65 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
66 cm1 * phi(i, j - 1, k, 2)) *
67 idx[1];
68
69 cp1 = Stencil::c20;
70 c = Stencil::c21;
71 cm1 = Stencil::c22;
72
73 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
74 cm1 * phi(i, j, k - 1, 0)) *
75 idx[2];
76 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
77 cm1 * phi(i, j, k - 1, 1)) *
78 idx[2];
79 wz = (cp1 * phi(i, j, k + 1, 2) + c * phi(i, j, k, 2) +
80 cm1 * phi(i, j, k - 1, 2)) *
81 idx[2];
82
83 strphi(i, j, k) = std::sqrt(
84 2.0_rt * std::pow(ux, 2.0_rt) +
85 2.0_rt * std::pow(vy, 2.0_rt) +
86 2.0_rt * std::pow(wz, 2.0_rt) + std::pow(uy + vx, 2.0_rt) +
87 std::pow(vz + wy, 2.0_rt) + std::pow(wx + uz, 2.0_rt));
88 });
89 }
90
91 FTypeOut& m_strphi;
92 const FTypeIn& m_phi;
93};
94
101template <typename FTypeIn, typename FTypeOut>
102inline void strainrate(FTypeOut& strphi, const FTypeIn& phi)
103{
104 BL_PROFILE("amr-wind::fvm::strainrate");
105 StrainRate<FTypeIn, FTypeOut> str(strphi, phi);
106 impl::apply(str, phi);
107}
108
114template <typename FType>
115inline std::unique_ptr<ScratchField> strainrate(const FType& phi)
116{
117 const std::string gname = phi.name() + "_strainrate";
118 auto strphi = phi.repo().create_scratch_field(gname, 1);
119 strainrate(*strphi, phi);
120 return strphi;
121}
122
123} // namespace amr_wind::fvm
124
125#endif /* STRAINRATE_H */
void strainrate(FTypeOut &strphi, const FTypeIn &phi)
Definition strainrate.H:102
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition strainrate.H:16
const FTypeIn & m_phi
Definition strainrate.H:92
void apply(const int lev, const amrex::MFIter &mfi) const
Definition strainrate.H:25
FTypeOut & m_strphi
Definition strainrate.H:91
StrainRate(FTypeOut &strphi, const FTypeIn &phi)
Definition strainrate.H:17