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