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

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/fvm/vorticity.H Source File
Kynema-SGF API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
vorticity.H
Go to the documentation of this file.
1#ifndef VORTICITY_H
2#define VORTICITY_H
3
4#include "src/fvm/fvm_utils.H"
5
6namespace kynema_sgf::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
14 Vorticity(FTypeOut& vortphi, const FTypeIn& phi)
15 : m_vort(vortphi), m_phi(phi)
16 {
17 AMREX_ALWAYS_ASSERT(AMREX_SPACEDIM == m_phi.num_comp());
18 AMREX_ALWAYS_ASSERT(m_phi.num_grow() > amrex::IntVect(0));
19 }
20
21 template <typename Stencil>
22 void apply(const int lev) const
23 {
24 const auto& geom = m_phi.repo().mesh().Geom(lev);
25 const auto& idx = geom.InvCellSizeArray();
26 const auto dlo = geom.Domain().smallEnd();
27 const auto dhi = geom.Domain().bigEnd();
28
29 auto& out_mf = m_vort(lev);
30 auto const& in_arrs = m_phi(lev).const_arrays();
31 auto const& out_arrs = out_mf.arrays();
32
33 amrex::ParallelFor(
34 out_mf, amrex::IntVect(0),
35 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
36 if (!Stencil::applies_to(i, j, k, dlo, dhi)) {
37 return;
38 }
39 auto const& phi = in_arrs[nbx];
40 auto const& vort = out_arrs[nbx];
41 const amrex::Real vx = (Stencil::c00 * phi(i + 1, j, k, 1) +
42 Stencil::c01 * phi(i, j, k, 1) +
43 Stencil::c02 * phi(i - 1, j, k, 1)) *
44 idx[0];
45 const amrex::Real wx = (Stencil::c00 * phi(i + 1, j, k, 2) +
46 Stencil::c01 * phi(i, j, k, 2) +
47 Stencil::c02 * phi(i - 1, j, k, 2)) *
48 idx[0];
49 const amrex::Real uy = (Stencil::c10 * phi(i, j + 1, k, 0) +
50 Stencil::c11 * phi(i, j, k, 0) +
51 Stencil::c12 * phi(i, j - 1, k, 0)) *
52 idx[1];
53 const amrex::Real wy = (Stencil::c10 * phi(i, j + 1, k, 2) +
54 Stencil::c11 * phi(i, j, k, 2) +
55 Stencil::c12 * phi(i, j - 1, k, 2)) *
56 idx[1];
57 const amrex::Real uz = (Stencil::c20 * phi(i, j, k + 1, 0) +
58 Stencil::c21 * phi(i, j, k, 0) +
59 Stencil::c22 * phi(i, j, k - 1, 0)) *
60 idx[2];
61 const amrex::Real vz = (Stencil::c20 * phi(i, j, k + 1, 1) +
62 Stencil::c21 * phi(i, j, k, 1) +
63 Stencil::c22 * phi(i, j, k - 1, 1)) *
64 idx[2];
65 vort(i, j, k, 0) = wy - vz;
66 vort(i, j, k, 1) = uz - wx;
67 vort(i, j, k, 2) = vx - uy;
68 });
69 }
70
71 FTypeOut& m_vort;
72 const FTypeIn& m_phi;
73};
74
81template <typename FTypeIn, typename FTypeOut>
82void vorticity(FTypeOut& vortphi, const FTypeIn& phi)
83{
84 BL_PROFILE("kynema-sgf::fvm::vorticity");
85 Vorticity<FTypeIn, FTypeOut> vort(vortphi, phi);
86 impl::apply(vort, phi);
87}
88
94template <typename FType>
95std::unique_ptr<ScratchField> vorticity(const FType& phi)
96{
97 const std::string gname = phi.name() + "_vorticity";
98 auto vortphi = phi.repo().create_scratch_field(gname, AMREX_SPACEDIM);
99 vorticity(*vortphi, phi);
100 return vortphi;
101}
102
103} // namespace kynema_sgf::fvm
104
105#endif /* VORTICITY_H */
void vorticity(FTypeOut &vortphi, const FTypeIn &phi)
Definition vorticity.H:82
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:20
Definition SchemeTraits.H:6
Definition vorticity.H:13
void apply(const int lev) const
Definition vorticity.H:22
FTypeOut & m_vort
Definition vorticity.H:71
const FTypeIn & m_phi
Definition vorticity.H:72
Vorticity(FTypeOut &vortphi, const FTypeIn &phi)
Definition vorticity.H:14