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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/fvm/vorticity.H Source File
AMR-Wind 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
5
6namespace amr_wind::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 amrex::MFIter& mfi) const
23 {
24 const auto& geom = m_phi.repo().mesh().Geom(lev);
25 const auto& idx = geom.InvCellSizeArray();
26 const auto& vort = m_vort(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(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
36 amrex::Real cp1, c, cm1, uy, uz, vx, vz, wx, wy;
37 cp1 = Stencil::c00;
38 c = Stencil::c01;
39 cm1 = Stencil::c02;
40
41 vx = (cp1 * phi(i + 1, j, k, 1) + c * phi(i, j, k, 1) +
42 cm1 * phi(i - 1, j, k, 1)) *
43 idx[0];
44 wx = (cp1 * phi(i + 1, j, k, 2) + c * phi(i, j, k, 2) +
45 cm1 * phi(i - 1, j, k, 2)) *
46 idx[0];
47
48 cp1 = Stencil::c10;
49 c = Stencil::c11;
50 cm1 = Stencil::c12;
51
52 uy = (cp1 * phi(i, j + 1, k, 0) + c * phi(i, j, k, 0) +
53 cm1 * phi(i, j - 1, k, 0)) *
54 idx[1];
55 wy = (cp1 * phi(i, j + 1, k, 2) + c * phi(i, j, k, 2) +
56 cm1 * phi(i, j - 1, k, 2)) *
57 idx[1];
58
59 cp1 = Stencil::c20;
60 c = Stencil::c21;
61 cm1 = Stencil::c22;
62
63 uz = (cp1 * phi(i, j, k + 1, 0) + c * phi(i, j, k, 0) +
64 cm1 * phi(i, j, k - 1, 0)) *
65 idx[2];
66 vz = (cp1 * phi(i, j, k + 1, 1) + c * phi(i, j, k, 1) +
67 cm1 * phi(i, j, k - 1, 1)) *
68 idx[2];
69
70 vort(i, j, k, 0) = wy - vz;
71 vort(i, j, k, 1) = uz - wx;
72 vort(i, j, k, 2) = vx - uy;
73 });
74 }
75
76 FTypeOut& m_vort;
77 const FTypeIn& m_phi;
78};
79
86template <typename FTypeIn, typename FTypeOut>
87void vorticity(FTypeOut& vortphi, const FTypeIn& phi)
88{
89 BL_PROFILE("amr-wind::fvm::vorticity");
90 Vorticity<FTypeIn, FTypeOut> vort(vortphi, phi);
91 impl::apply(vort, phi);
92}
93
99template <typename FType>
100std::unique_ptr<ScratchField> vorticity(const FType& phi)
101{
102 const std::string gname = phi.name() + "_vorticity";
103 auto vortphi = phi.repo().create_scratch_field(gname, AMREX_SPACEDIM);
104 vorticity(*vortphi, phi);
105 return vortphi;
106}
107
108} // namespace amr_wind::fvm
109
110#endif /* VORTICITY_H */
void vorticity(FTypeOut &vortphi, const FTypeIn &phi)
Definition vorticity.H:87
AMREX_INLINE void apply(const FvmOp &fvmop, const FType &fld)
Definition fvm_utils.H:14
Definition SchemeTraits.H:6
Definition vorticity.H:13
const FTypeIn & m_phi
Definition vorticity.H:77
Vorticity(FTypeOut &vortphi, const FTypeIn &phi)
Definition vorticity.H:14
void apply(const int lev, const amrex::MFIter &mfi) const
Definition vorticity.H:22
FTypeOut & m_vort
Definition vorticity.H:76