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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/vorticity_mag.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
vorticity_mag.H
Go to the documentation of this file.
1#ifndef VORTICITYMAG_H
2#define VORTICITYMAG_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 VorticityMag(FTypeOut& vortmagphi, const FTypeIn& phi)
18 : m_vortmagphi(vortmagphi), 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& vortmagphi = m_vortmagphi(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(
39 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
40 amrex::Real cp1, c, cm1, uy, uz, vx, vz, wx, wy;
41 cp1 = Stencil::c00;
42 c = Stencil::c01;
43 cm1 = Stencil::c02;
44
45 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
46 cm1 * phi(i - 1, j, k, 1)) *
47 idx[0];
48 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
49 cm1 * phi(i - 1, j, k, 2)) *
50 idx[0];
51
52 cp1 = Stencil::c10;
53 c = Stencil::c11;
54 cm1 = Stencil::c12;
55
56 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
57 cm1 * phi(i, j - 1, k, 0)) *
58 idx[1];
59 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
60 cm1 * phi(i, j - 1, k, 2)) *
61 idx[1];
62
63 cp1 = Stencil::c20;
64 c = Stencil::c21;
65 cm1 = Stencil::c22;
66
67 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
68 cm1 * phi(i, j, k - 1, 0)) *
69 idx[2];
70 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
71 cm1 * phi(i, j, k - 1, 1)) *
72 idx[2];
73
74 vortmagphi(i, j, k) = std::sqrt(
75 std::pow(uy - vx, 2.0_rt) + std::pow(vz - wy, 2.0_rt) +
76 std::pow(wx - uz, 2.0_rt));
77 });
78 }
79
80 FTypeOut& m_vortmagphi;
81 const FTypeIn& m_phi;
82};
83
90template <typename FTypeIn, typename FTypeOut>
91inline void vorticity_mag(FTypeOut& vortmagphi, const FTypeIn& phi)
92{
93 BL_PROFILE("amr-wind::fvm::vorticity_mag");
94 VorticityMag<FTypeIn, FTypeOut> vortmag(vortmagphi, phi);
95 impl::apply(vortmag, phi);
96}
97
103template <typename FType>
104inline std::unique_ptr<ScratchField> vorticity_mag(const FType& phi)
105{
106 const std::string gname = phi.name() + "_vorticity_mag";
107 auto vortmagphi = phi.repo().create_scratch_field(gname, 1);
108 vorticity_mag(*vortmagphi, phi);
109 return vortmagphi;
110}
111
112} // namespace amr_wind::fvm
113
114#endif /* VORTICITYMAG_H */
void vorticity_mag(FTypeOut &vortmagphi, const FTypeIn &phi)
Definition vorticity_mag.H:91
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition vorticity_mag.H:16
void apply(const int lev, const amrex::MFIter &mfi) const
Definition vorticity_mag.H:25
const FTypeIn & m_phi
Definition vorticity_mag.H:81
FTypeOut & m_vortmagphi
Definition vorticity_mag.H:80
VorticityMag(FTypeOut &vortmagphi, const FTypeIn &phi)
Definition vorticity_mag.H:17