/home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/CompRHSOps.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/CompRHSOps.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
CompRHSOps.H
Go to the documentation of this file.
1#ifndef COMPRHSOPS_H
2#define COMPRHSOPS_H
3
7
8namespace amr_wind::pde {
9
16template <typename PDE, typename Scheme>
18{
19 explicit ComputeRHSOp(PDEFields& fields_in)
20 : fields(fields_in), density(fields_in.repo.get_field("density"))
21 {}
22
30 const DiffusionType difftype, const amrex::Real dt, bool mesh_mapping)
31 {
32 amrex::Real factor = 0.0;
33 switch (difftype) {
35 factor = 1.0;
36 break;
37
39 factor = 0.5;
40 break;
41
43 factor = 0.0;
44 break;
45
46 default:
47 amrex::Abort("Invalid diffusion type");
48 }
49
50 // Field states for diffusion and advection terms. In Godunov scheme
51 // these terms only have one state.
52 auto fstate = std::is_same<Scheme, fvm::Godunov>::value
55
56 const int nlevels = fields.repo.num_active_levels();
57
58 // for RHS evaluation velocity field should be in stretched space
59 auto& field = fields.field;
60 if (field.in_uniform_space() && mesh_mapping) {
61 field.to_stretched_space();
62 }
63 auto& field_old = field.state(FieldState::Old);
64 if (field_old.in_uniform_space() && mesh_mapping) {
65 field_old.to_stretched_space();
66 }
67
68 auto& den_new = density.state(FieldState::New);
69 auto& den_old = density.state(FieldState::Old);
70 auto& src_term = fields.src_term;
71 auto& diff_term = fields.diff_term.state(fstate);
72 auto& conv_term = fields.conv_term.state(fstate);
73 auto& mask_cell = fields.repo.get_int_field("mask_cell");
74 Field const* mesh_detJ =
76 : nullptr;
77
78 for (int lev = 0; lev < nlevels; ++lev) {
79 const auto& fld_arrs = field(lev).arrays();
80 const auto& fld_o_arrs = field_old(lev).const_arrays();
81 const auto& rho_o_arrs = den_old(lev).const_arrays();
82 const auto& rho_arrs = den_new(lev).const_arrays();
83 const auto& src_arrs = src_term(lev).const_arrays();
84 const auto& diff_arrs = diff_term(lev).const_arrays();
85 const auto& ddt_o_arrs = conv_term(lev).const_arrays();
86 const auto& imask_arrs = mask_cell(lev).const_arrays();
87 const auto& detJ_arrs =
88 mesh_mapping ? ((*mesh_detJ)(lev).const_arrays())
89 : amrex::MultiArray4<amrex::Real const>();
90
91 if (PDE::multiply_rho) {
92 // Remove multiplication by density as it will be added back
93 // in solver
94 amrex::ParallelFor(
95 field(lev), amrex::IntVect(0), PDE::ndim,
96 [=] AMREX_GPU_DEVICE(
97 int nbx, int i, int j, int k, int n) noexcept {
98 amrex::Real det_j =
99 mesh_mapping ? (detJ_arrs[nbx](i, j, k)) : 1.0;
100
101 fld_arrs[nbx](i, j, k, n) =
102 rho_o_arrs[nbx](i, j, k) * det_j *
103 fld_o_arrs[nbx](i, j, k, n) +
104 static_cast<amrex::Real>(imask_arrs[nbx](i, j, k)) *
105 dt *
106 (ddt_o_arrs[nbx](i, j, k, n) +
107 det_j * src_arrs[nbx](i, j, k, n) +
108 factor * diff_arrs[nbx](i, j, k, n));
109
110 fld_arrs[nbx](i, j, k, n) /= rho_arrs[nbx](i, j, k);
111
112 if (difftype == DiffusionType::Explicit) {
113 fld_arrs[nbx](i, j, k, n) /= det_j;
114 }
115 });
116 } else {
117 amrex::ParallelFor(
118 field(lev), amrex::IntVect(0), PDE::ndim,
119 [=] AMREX_GPU_DEVICE(
120 int nbx, int i, int j, int k, int n) noexcept {
121 amrex::Real det_j =
122 mesh_mapping ? (detJ_arrs[nbx](i, j, k)) : 1.0;
123
124 fld_arrs[nbx](i, j, k, n) =
125 det_j * fld_o_arrs[nbx](i, j, k, n) +
126 static_cast<amrex::Real>(imask_arrs[nbx](i, j, k)) *
127 dt *
128 (ddt_o_arrs[nbx](i, j, k, n) +
129 det_j * src_arrs[nbx](i, j, k, n) +
130 factor * diff_arrs[nbx](i, j, k, n));
131
132 if (difftype == DiffusionType::Explicit) {
133 fld_arrs[nbx](i, j, k, n) /= det_j;
134 }
135 });
136 }
137 }
138 amrex::Gpu::streamSynchronize();
139 }
140
148 const DiffusionType difftype, const amrex::Real dt, bool mesh_mapping)
149 {
150 amrex::Real ofac = 0.0;
151 amrex::Real nfac = 0.0;
152 switch (difftype) {
154 ofac = 0.5;
155 nfac = 0.5;
156 break;
157
159 ofac = 0.5;
160 nfac = 0.0;
161 break;
162
164 ofac = 0.0;
165 nfac = 0.0;
166 break;
167
168 default:
169 amrex::Abort("Invalid diffusion type");
170 }
171
172 const int nlevels = fields.repo.num_active_levels();
173
174 // for RHS evaluation velocity field should be in stretched space
175 auto& field = fields.field;
176 if (field.in_uniform_space() && mesh_mapping) {
177 field.to_stretched_space();
178 }
179 auto& field_old = field.state(FieldState::Old);
180 if (field_old.in_uniform_space() && mesh_mapping) {
181 field_old.to_stretched_space();
182 }
183
184 auto& den_new = density.state(FieldState::New);
185 auto& den_old = density.state(FieldState::Old);
186 auto& src_term = fields.src_term;
187 auto& diff_term = fields.diff_term;
188 auto& conv_term = fields.conv_term;
189 auto& diff_term_old = fields.diff_term.state(FieldState::Old);
190 auto& conv_term_old = fields.conv_term.state(FieldState::Old);
191 auto& mask_cell = fields.repo.get_int_field("mask_cell");
192 Field const* mesh_detJ =
194 : nullptr;
195
196 for (int lev = 0; lev < nlevels; ++lev) {
197 const auto& fld_arrs = field(lev).arrays();
198 const auto& fld_o_arrs = field_old(lev).const_arrays();
199 const auto& rho_o_arrs = den_old(lev).const_arrays();
200 const auto& rho_arrs = den_new(lev).const_arrays();
201 const auto& src_arrs = src_term(lev).const_arrays();
202 const auto& diff_arrs = diff_term(lev).const_arrays();
203 const auto& ddt_arrs = conv_term(lev).const_arrays();
204 const auto& diff_o_arrs = diff_term_old(lev).const_arrays();
205 const auto& ddt_o_arrs = conv_term_old(lev).const_arrays();
206 const auto& imask_arrs = mask_cell(lev).const_arrays();
207 const auto& detJ_arrs =
208 mesh_mapping ? ((*mesh_detJ)(lev).const_arrays())
209 : amrex::MultiArray4<amrex::Real const>();
210
211 if (PDE::multiply_rho) {
212 // Remove multiplication by density as it will be added back
213 // in solver
214 amrex::ParallelFor(
215 field(lev), amrex::IntVect(0), PDE::ndim,
216 [=] AMREX_GPU_DEVICE(
217 int nbx, int i, int j, int k, int n) noexcept {
218 amrex::Real det_j =
219 mesh_mapping ? (detJ_arrs[nbx](i, j, k)) : 1.0;
220
221 fld_arrs[nbx](i, j, k, n) =
222 rho_o_arrs[nbx](i, j, k) * det_j *
223 fld_o_arrs[nbx](i, j, k, n) +
224 static_cast<amrex::Real>(imask_arrs[nbx](i, j, k)) *
225 dt *
226 (0.5 * (ddt_o_arrs[nbx](i, j, k, n) +
227 ddt_arrs[nbx](i, j, k, n)) +
228 ofac * diff_o_arrs[nbx](i, j, k, n) +
229 nfac * diff_arrs[nbx](i, j, k, n) +
230 det_j * src_arrs[nbx](i, j, k, n));
231
232 fld_arrs[nbx](i, j, k, n) /= rho_arrs[nbx](i, j, k);
233
234 if (difftype == DiffusionType::Explicit) {
235 fld_arrs[nbx](i, j, k, n) /= det_j;
236 }
237 });
238 } else {
239 amrex::ParallelFor(
240 field(lev), amrex::IntVect(0), PDE::ndim,
241 [=] AMREX_GPU_DEVICE(
242 int nbx, int i, int j, int k, int n) noexcept {
243 amrex::Real det_j =
244 mesh_mapping ? (detJ_arrs[nbx](i, j, k)) : 1.0;
245
246 fld_arrs[nbx](i, j, k, n) =
247 det_j * fld_o_arrs[nbx](i, j, k, n) +
248 static_cast<amrex::Real>(imask_arrs[nbx](i, j, k)) *
249 dt *
250 (0.5 * (ddt_o_arrs[nbx](i, j, k, n) +
251 ddt_arrs[nbx](i, j, k, n)) +
252 ofac * diff_o_arrs[nbx](i, j, k, n) +
253 nfac * diff_arrs[nbx](i, j, k, n) +
254 det_j * src_arrs[nbx](i, j, k, n));
255
256 if (difftype == DiffusionType::Explicit) {
257 fld_arrs[nbx](i, j, k, n) /= det_j;
258 }
259 });
260 }
261 }
262 amrex::Gpu::streamSynchronize();
263 }
264
265 // data members
268};
269
270} // namespace amr_wind::pde
271
272#endif /* COMPRHSOPS_H */
Definition Field.H:116
Field & state(const FieldState fstate)
Return field at a different time state.
Definition Field.cpp:114
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
Field & get_mesh_mapping_det_j(FieldLoc floc) const
Definition FieldRepo.cpp:190
IntField & get_int_field(const std::string &name, const FieldState fstate=FieldState::New) const
Return a reference to an integer field.
Definition FieldRepo.cpp:276
@ CELL
Cell-centered (default)
@ New
Same as FieldState::NP1.
@ Old
Same as FieldState::N.
DiffusionType
Definition incflo_enums.H:4
Definition AdvOp_Godunov.H:16
Definition CompRHSOps.H:18
ComputeRHSOp(PDEFields &fields_in)
Definition CompRHSOps.H:19
void predictor_rhs(const DiffusionType difftype, const amrex::Real dt, bool mesh_mapping)
Definition CompRHSOps.H:29
void corrector_rhs(const DiffusionType difftype, const amrex::Real dt, bool mesh_mapping)
Definition CompRHSOps.H:147
PDEFields & fields
Definition CompRHSOps.H:266
Field & density
Definition CompRHSOps.H:267
Definition PDEFields.H:27
Field & src_term
Source term for this PDE.
Definition PDEFields.H:39
Field & diff_term
Diffusion term for this PDE.
Definition PDEFields.H:41
Field & conv_term
Convective term for this PDE.
Definition PDEFields.H:43
FieldRepo & repo
Reference to the field repository instance.
Definition PDEFields.H:31
Field & field
Solution variable (e.g., velocity, temperature)
Definition PDEFields.H:34