/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
6namespace amr_wind::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
15 FTypeOut& qcritphi, const FTypeIn& phi, const bool nondim = false)
16 : m_qcritphi(qcritphi), m_phi(phi), m_nondim(nondim)
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& qcritphi = m_qcritphi(lev).array(mfi);
28 const auto& phi = m_phi(lev).const_array(mfi);
29 const bool nondim = m_nondim;
30
31 const auto& bx_in = mfi.tilebox();
32 const auto& bx = Stencil::box(bx_in, geom);
33 if (bx.isEmpty()) {
34 return;
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 const amrex::Real S2 =
82 std::pow(ux, 2) + std::pow(vy, 2) + std::pow(wz, 2) +
83 0.5 * std::pow(uy + vx, 2) + 0.5 * std::pow(vz + wy, 2) +
84 0.5 * std::pow(wx + uz, 2);
85
86 const amrex::Real W2 = 0.5 * std::pow(uy - vx, 2) +
87 0.5 * std::pow(vz - wy, 2) +
88 0.5 * std::pow(wx - uz, 2);
89 if (nondim) {
90 qcritphi(i, j, k) =
91 0.5 * (W2 / amrex::max(1e-14, S2) - 1.0);
92 } else {
93 qcritphi(i, j, k) = 0.5 * (W2 - S2);
94 }
95 });
96 }
97
98 FTypeOut& m_qcritphi;
99 const FTypeIn& m_phi;
101};
102
110template <typename FTypeIn, typename FTypeOut>
111inline void
112q_criterion(FTypeOut& qcritphi, const FTypeIn& phi, const bool nondim = false)
113{
114 BL_PROFILE("amr-wind::fvm::q_criterion");
115 Qcriterion<FTypeIn, FTypeOut> qcriterion(qcritphi, phi, nondim);
116 impl::apply(qcriterion, phi);
117}
118
124template <typename FType>
125inline std::unique_ptr<ScratchField> q_criterion(const FType& phi)
126{
127 const std::string gname = phi.name() + "_q_criterion";
128 auto qcritphi = phi.repo().create_scratch_field(gname, 1);
129 q_criterion(*qcritphi, phi);
130 return qcritphi;
131}
132
133} // namespace amr_wind::fvm
134
135#endif /* QCRITERION_H */
void q_criterion(FTypeOut &qcritphi, const FTypeIn &phi, const bool nondim=false)
Definition qcriterion.H:112
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition qcriterion.H:13
FTypeOut & m_qcritphi
Definition qcriterion.H:98
bool m_nondim
Definition qcriterion.H:100
const FTypeIn & m_phi
Definition qcriterion.H:99
Qcriterion(FTypeOut &qcritphi, const FTypeIn &phi, const bool nondim=false)
Definition qcriterion.H:14
void apply(const int lev, const amrex::MFIter &mfi) const
Definition qcriterion.H:23