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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/qcriterion.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
qcriterion.H
Go to the documentation of this file.
1#ifndef QCRITERION_H
2#define QCRITERION_H
3
5#include "AMReX_REAL.H"
6
7using namespace amrex::literals;
8
9namespace amr_wind::fvm {
10
14template <typename FTypeIn, typename FTypeOut>
16{
18 FTypeOut& qcritphi, const FTypeIn& phi, const bool nondim = false)
19 : m_qcritphi(qcritphi), m_phi(phi), m_nondim(nondim)
20 {
21 AMREX_ALWAYS_ASSERT(AMREX_SPACEDIM == m_phi.num_comp());
22 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
23 }
24
25 template <typename Stencil>
26 void apply(const int lev, const amrex::MFIter& mfi) const
27 {
28 const auto& geom = m_phi.repo().mesh().Geom(lev);
29 const auto& idx = geom.InvCellSizeArray();
30 const auto& qcritphi = m_qcritphi(lev).array(mfi);
31 const auto& phi = m_phi(lev).const_array(mfi);
32 const bool nondim = m_nondim;
33
34 const auto& bx_in = mfi.tilebox();
35 const auto& bx = Stencil::box(bx_in, geom);
36 if (bx.isEmpty()) {
37 return;
38 }
39 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
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 const amrex::Real S2 = std::pow(ux, 2.0_rt) + std::pow(vy, 2.0_rt) +
84 std::pow(wz, 2.0_rt) +
85 (0.5_rt * std::pow(uy + vx, 2.0_rt)) +
86 (0.5_rt * std::pow(vz + wy, 2.0_rt)) +
87 (0.5_rt * std::pow(wx + uz, 2.0_rt));
88
89 const amrex::Real W2 = (0.5_rt * std::pow(uy - vx, 2.0_rt)) +
90 (0.5_rt * std::pow(vz - wy, 2.0_rt)) +
91 (0.5_rt * std::pow(wx - uz, 2.0_rt));
92 if (nondim) {
93 qcritphi(i, j, k) =
94 0.5_rt *
95 (W2 / amrex::max<amrex::Real>(
96 std::numeric_limits<amrex::Real>::epsilon() *
97 1.0e2_rt,
98 S2) -
99 1.0_rt);
100 } else {
101 qcritphi(i, j, k) = 0.5_rt * (W2 - S2);
102 }
103 });
104 }
105
106 FTypeOut& m_qcritphi;
107 const FTypeIn& m_phi;
109};
110
118template <typename FTypeIn, typename FTypeOut>
120 FTypeOut& qcritphi, const FTypeIn& phi, const bool nondim = false)
121{
122 BL_PROFILE("amr-wind::fvm::q_criterion");
123 Qcriterion<FTypeIn, FTypeOut> qcriterion(qcritphi, phi, nondim);
124 impl::apply(qcriterion, phi);
125}
126
132template <typename FType>
133std::unique_ptr<ScratchField> q_criterion(const FType& phi)
134{
135 const std::string gname = phi.name() + "_q_criterion";
136 auto qcritphi = phi.repo().create_scratch_field(gname, 1);
137 q_criterion(*qcritphi, phi);
138 return qcritphi;
139}
140
141} // namespace amr_wind::fvm
142
143#endif /* QCRITERION_H */
void q_criterion(FTypeOut &qcritphi, const FTypeIn &phi, const bool nondim=false)
Definition qcriterion.H:119
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition qcriterion.H:16
FTypeOut & m_qcritphi
Definition qcriterion.H:106
bool m_nondim
Definition qcriterion.H:108
const FTypeIn & m_phi
Definition qcriterion.H:107
Qcriterion(FTypeOut &qcritphi, const FTypeIn &phi, const bool nondim=false)
Definition qcriterion.H:17
void apply(const int lev, const amrex::MFIter &mfi) const
Definition qcriterion.H:26