/home/runner/work/kynema-sgf/kynema-sgf/src/fvm/qcriterion.H Source File

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/fvm/qcriterion.H Source File
Kynema-SGF 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
4#include "src/fvm/fvm_utils.H"
5#include "AMReX_REAL.H"
6
7using namespace amrex::literals;
8
9namespace kynema_sgf::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
27 {
28 const auto& geom = m_phi.repo().mesh().Geom(lev);
29 const auto& idx = geom.InvCellSizeArray();
30 const auto dlo = geom.Domain().smallEnd();
31 const auto dhi = geom.Domain().bigEnd();
32 const bool nondim = m_nondim;
33
34 auto& out_mf = m_qcritphi(lev);
35 auto const& in_arrs = m_phi(lev).const_arrays();
36 auto const& out_arrs = out_mf.arrays();
37
38 amrex::ParallelFor(
39 out_mf, amrex::IntVect(0),
40 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
41 if (!Stencil::applies_to(i, j, k, dlo, dhi)) {
42 return;
43 }
44 auto const& phi = in_arrs[nbx];
45 auto const& qcritphi = out_arrs[nbx];
46 const amrex::Real ux = (Stencil::c00 * phi(i + 1, j, k, 0) +
47 Stencil::c01 * phi(i, j, k, 0) +
48 Stencil::c02 * phi(i - 1, j, k, 0)) *
49 idx[0];
50 const amrex::Real vx = (Stencil::c00 * phi(i + 1, j, k, 1) +
51 Stencil::c01 * phi(i, j, k, 1) +
52 Stencil::c02 * phi(i - 1, j, k, 1)) *
53 idx[0];
54 const amrex::Real wx = (Stencil::c00 * phi(i + 1, j, k, 2) +
55 Stencil::c01 * phi(i, j, k, 2) +
56 Stencil::c02 * phi(i - 1, j, k, 2)) *
57 idx[0];
58 const amrex::Real uy = (Stencil::c10 * phi(i, j + 1, k, 0) +
59 Stencil::c11 * phi(i, j, k, 0) +
60 Stencil::c12 * phi(i, j - 1, k, 0)) *
61 idx[1];
62 const amrex::Real vy = (Stencil::c10 * phi(i, j + 1, k, 1) +
63 Stencil::c11 * phi(i, j, k, 1) +
64 Stencil::c12 * phi(i, j - 1, k, 1)) *
65 idx[1];
66 const amrex::Real wy = (Stencil::c10 * phi(i, j + 1, k, 2) +
67 Stencil::c11 * phi(i, j, k, 2) +
68 Stencil::c12 * phi(i, j - 1, k, 2)) *
69 idx[1];
70 const amrex::Real uz = (Stencil::c20 * phi(i, j, k + 1, 0) +
71 Stencil::c21 * phi(i, j, k, 0) +
72 Stencil::c22 * phi(i, j, k - 1, 0)) *
73 idx[2];
74 const amrex::Real vz = (Stencil::c20 * phi(i, j, k + 1, 1) +
75 Stencil::c21 * phi(i, j, k, 1) +
76 Stencil::c22 * phi(i, j, k - 1, 1)) *
77 idx[2];
78 const amrex::Real wz = (Stencil::c20 * phi(i, j, k + 1, 2) +
79 Stencil::c21 * phi(i, j, k, 2) +
80 Stencil::c22 * phi(i, j, k - 1, 2)) *
81 idx[2];
82 const amrex::Real S2 = utils::powi(ux, 2) + utils::powi(vy, 2) +
83 utils::powi(wz, 2) +
84 (0.5_rt * utils::powi(uy + vx, 2)) +
85 (0.5_rt * utils::powi(vz + wy, 2)) +
86 (0.5_rt * utils::powi(wx + uz, 2));
87 const amrex::Real W2 = (0.5_rt * utils::powi(uy - vx, 2)) +
88 (0.5_rt * utils::powi(vz - wy, 2)) +
89 (0.5_rt * utils::powi(wx - uz, 2));
90 if (nondim) {
91 qcritphi(i, j, k) =
92 0.5_rt *
93 (W2 / amrex::max<amrex::Real>(
94 std::numeric_limits<amrex::Real>::epsilon() *
95 1.0e2_rt,
96 S2) -
97 1.0_rt);
98 } else {
99 qcritphi(i, j, k) = 0.5_rt * (W2 - S2);
100 }
101 });
102 }
103
104 FTypeOut& m_qcritphi;
105 const FTypeIn& m_phi;
107};
108
116template <typename FTypeIn, typename FTypeOut>
118 FTypeOut& qcritphi, const FTypeIn& phi, const bool nondim = false)
119{
120 BL_PROFILE("kynema-sgf::fvm::q_criterion");
121 Qcriterion<FTypeIn, FTypeOut> qcriterion(qcritphi, phi, nondim);
122 impl::apply(qcriterion, phi);
123}
124
130template <typename FType>
131std::unique_ptr<ScratchField> q_criterion(const FType& phi)
132{
133 const std::string gname = phi.name() + "_q_criterion";
134 auto qcritphi = phi.repo().create_scratch_field(gname, 1);
135 q_criterion(*qcritphi, phi);
136 return qcritphi;
137}
138
139} // namespace kynema_sgf::fvm
140
141#endif /* QCRITERION_H */
void q_criterion(FTypeOut &qcritphi, const FTypeIn &phi, const bool nondim=false)
Definition qcriterion.H:117
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:20
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 qcriterion.H:16
FTypeOut & m_qcritphi
Definition qcriterion.H:104
void apply(const int lev) const
Definition qcriterion.H:26
bool m_nondim
Definition qcriterion.H:106
const FTypeIn & m_phi
Definition qcriterion.H:105
Qcriterion(FTypeOut &qcritphi, const FTypeIn &phi, const bool nondim=false)
Definition qcriterion.H:17