/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
6namespace amr_wind::fvm {
7
11template <typename FTypeIn, typename FTypeOut>
13{
14 VorticityMag(FTypeOut& vortmagphi, const FTypeIn& phi)
15 : m_vortmagphi(vortmagphi), 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 amrex::MFIter& mfi) const
23 {
24 const auto& geom = m_phi.repo().mesh().Geom(lev);
25 const auto& idx = geom.InvCellSizeArray();
26 const auto& vortmagphi = m_vortmagphi(lev).array(mfi);
27 const auto& phi = m_phi(lev).const_array(mfi);
28
29 const auto& bx_in = mfi.tilebox();
30 const auto& bx = Stencil::box(bx_in, geom);
31 if (bx.isEmpty()) {
32 return;
33 }
34
35 amrex::ParallelFor(
36 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
37 amrex::Real cp1, c, cm1, uy, uz, vx, vz, wx, wy;
38 cp1 = Stencil::c00;
39 c = Stencil::c01;
40 cm1 = Stencil::c02;
41
42 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
43 cm1 * phi(i - 1, j, k, 1)) *
44 idx[0];
45 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
46 cm1 * phi(i - 1, j, k, 2)) *
47 idx[0];
48
49 cp1 = Stencil::c10;
50 c = Stencil::c11;
51 cm1 = Stencil::c12;
52
53 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
54 cm1 * phi(i, j - 1, k, 0)) *
55 idx[1];
56 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
57 cm1 * phi(i, j - 1, k, 2)) *
58 idx[1];
59
60 cp1 = Stencil::c20;
61 c = Stencil::c21;
62 cm1 = Stencil::c22;
63
64 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
65 cm1 * phi(i, j, k - 1, 0)) *
66 idx[2];
67 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
68 cm1 * phi(i, j, k - 1, 1)) *
69 idx[2];
70
71 vortmagphi(i, j, k) = sqrt(
72 std::pow(uy - vx, 2) + std::pow(vz - wy, 2) +
73 std::pow(wx - uz, 2));
74 });
75 }
76
77 FTypeOut& m_vortmagphi;
78 const FTypeIn& m_phi;
79};
80
87template <typename FTypeIn, typename FTypeOut>
88inline void vorticity_mag(FTypeOut& vortmagphi, const FTypeIn& phi)
89{
90 BL_PROFILE("amr-wind::fvm::vorticity_mag");
91 VorticityMag<FTypeIn, FTypeOut> vortmag(vortmagphi, phi);
92 impl::apply(vortmag, phi);
93}
94
100template <typename FType>
101inline std::unique_ptr<ScratchField> vorticity_mag(const FType& phi)
102{
103 const std::string gname = phi.name() + "_vorticity_mag";
104 auto vortmagphi = phi.repo().create_scratch_field(gname, 1);
105 vorticity_mag(*vortmagphi, phi);
106 return vortmagphi;
107}
108
109} // namespace amr_wind::fvm
110
111#endif /* VORTICITYMAG_H */
void vorticity_mag(FTypeOut &vortmagphi, const FTypeIn &phi)
Definition vorticity_mag.H:88
void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition vorticity_mag.H:13
void apply(const int lev, const amrex::MFIter &mfi) const
Definition vorticity_mag.H:22
const FTypeIn & m_phi
Definition vorticity_mag.H:78
FTypeOut & m_vortmagphi
Definition vorticity_mag.H:77
VorticityMag(FTypeOut &vortmagphi, const FTypeIn &phi)
Definition vorticity_mag.H:14