/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(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
39 amrex::Real cp1, c, cm1, uy, uz, vx, vz, wx, wy;
40 cp1 = Stencil::c00;
41 c = Stencil::c01;
42 cm1 = Stencil::c02;
43
44 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
45 cm1 * phi(i - 1, j, k, 1)) *
46 idx[0];
47 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
48 cm1 * phi(i - 1, j, k, 2)) *
49 idx[0];
50
51 cp1 = Stencil::c10;
52 c = Stencil::c11;
53 cm1 = Stencil::c12;
54
55 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
56 cm1 * phi(i, j - 1, k, 0)) *
57 idx[1];
58 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
59 cm1 * phi(i, j - 1, k, 2)) *
60 idx[1];
61
62 cp1 = Stencil::c20;
63 c = Stencil::c21;
64 cm1 = Stencil::c22;
65
66 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
67 cm1 * phi(i, j, k - 1, 0)) *
68 idx[2];
69 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
70 cm1 * phi(i, j, k - 1, 1)) *
71 idx[2];
72
73 vortmagphi(i, j, k) = std::sqrt(
74 std::pow(uy - vx, 2.0_rt) + std::pow(vz - wy, 2.0_rt) +
75 std::pow(wx - uz, 2.0_rt));
76 });
77 }
78
79 FTypeOut& m_vortmagphi;
80 const FTypeIn& m_phi;
81};
82
89template <typename FTypeIn, typename FTypeOut>
90void vorticity_mag(FTypeOut& vortmagphi, const FTypeIn& phi)
91{
92 BL_PROFILE("amr-wind::fvm::vorticity_mag");
93 VorticityMag<FTypeIn, FTypeOut> vortmag(vortmagphi, phi);
94 impl::apply(vortmag, phi);
95}
96
102template <typename FType>
103std::unique_ptr<ScratchField> vorticity_mag(const FType& phi)
104{
105 const std::string gname = phi.name() + "_vorticity_mag";
106 auto vortmagphi = phi.repo().create_scratch_field(gname, 1);
107 vorticity_mag(*vortmagphi, phi);
108 return vortmagphi;
109}
110
111} // namespace amr_wind::fvm
112
113#endif /* VORTICITYMAG_H */
void vorticity_mag(FTypeOut &vortmagphi, const FTypeIn &phi)
Definition vorticity_mag.H:90
AMREX_INLINE 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:80
FTypeOut & m_vortmagphi
Definition vorticity_mag.H:79
VorticityMag(FTypeOut &vortmagphi, const FTypeIn &phi)
Definition vorticity_mag.H:17