/home/runner/work/amr-wind/amr-wind/amr-wind/core/Field.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/core/Field.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
Field.H
Go to the documentation of this file.
1#ifndef FIELD_H
2#define FIELD_H
3
4#include <string>
5#include <memory>
6#include <unordered_map>
7#include <algorithm>
8
9#include "AMReX_MultiFab.H"
10#include "AMReX_BCRec.H"
11#include "AMReX_Gpu.H"
12
16
17namespace amr_wind {
18
19class Field;
20class FieldRepo;
22class FieldBCIface;
23class SimTime;
24
29{
31 static constexpr int max_field_states = 5;
32
34 std::string basename, int ncomp, int ngrow, int nstates, FieldLoc floc);
35
37
39 bool bc_initialized();
40
42 void copy_bc_to_device();
43
45 std::string m_basename;
46
49
51 amrex::IntVect m_ngrow;
52
55
58
61 amrex::GpuArray<BC, AMREX_SPACEDIM * 2> m_bc_type;
62 amrex::Vector<amrex::Vector<amrex::Real>> m_bc_values;
63 amrex::GpuArray<const amrex::Real*, AMREX_SPACEDIM * 2> m_bc_values_d;
64 amrex::Gpu::DeviceVector<amrex::Real> m_bc_values_dview;
65
66 amrex::Vector<amrex::BCRec> m_bcrec;
67 amrex::Gpu::DeviceVector<amrex::BCRec> m_bcrec_d;
69
71 amrex::Vector<Field*> m_states;
72
74 std::unique_ptr<FieldFillPatchOpsBase> m_fillpatch_op;
75
77 amrex::Vector<std::unique_ptr<FieldBCIface>> m_bc_func;
78
81};
82
111class Field
112{
113public:
114 friend class FieldRepo;
115
116 Field(const Field&) = delete;
117 Field& operator=(const Field&) = delete;
119
121 [[nodiscard]] const std::string& name() const { return m_name; }
122
124 [[nodiscard]] const std::string& base_name() const
125 {
126 return m_info->m_basename;
127 }
128
130 [[nodiscard]] unsigned id() const { return m_id; }
131
133 [[nodiscard]] int num_comp() const { return m_info->m_ncomp; }
134
136 [[nodiscard]] const amrex::IntVect& num_grow() const
137 {
138 return m_info->m_ngrow;
139 }
140
142 [[nodiscard]] int num_time_states() const { return m_info->m_nstates; }
143
145 [[nodiscard]] int num_states() const
146 {
147 int ns = static_cast<int>(std::count_if(
148 m_info->m_states.begin(), m_info->m_states.end(),
149 [](const Field* ff) { return (ff != nullptr); }));
150 AMREX_ASSERT(num_time_states() <= ns);
151 return ns;
152 }
153
155 [[nodiscard]] FieldLoc field_location() const { return m_info->m_floc; }
156
158 [[nodiscard]] FieldState field_state() const { return m_state; }
159
161 [[nodiscard]] FieldRepo& repo() const { return m_repo; }
162
164 [[nodiscard]] bool fillpatch_on_regrid() const
165 {
167 }
168
170 [[nodiscard]] bool query_state(const FieldState fstate) const
171 {
172 const int fid = static_cast<int>(fstate);
173 return (m_info->m_states[fid] != nullptr);
174 }
175
176 [[nodiscard]] const amrex::GpuArray<BC, AMREX_SPACEDIM * 2>& bc_type() const
177 {
178 return m_info->m_bc_type;
179 }
180
181 amrex::GpuArray<BC, AMREX_SPACEDIM * 2>& bc_type()
182 {
183 return m_info->m_bc_type;
184 }
185
187 amrex::Vector<amrex::Vector<amrex::Real>>& bc_values()
188 {
189 return m_info->m_bc_values;
190 }
191
193 [[nodiscard]] amrex::Vector<amrex::BCRec>& bcrec() const
194 {
195 return m_info->m_bcrec;
196 }
197
198 [[nodiscard]] const amrex::GpuArray<const amrex::Real*, AMREX_SPACEDIM * 2>&
200 {
201 return m_info->m_bc_values_d;
202 }
203
204 [[nodiscard]] const amrex::Gpu::DeviceVector<amrex::BCRec>&
206 {
207 return m_info->m_bcrec_d;
208 }
209
214 void copy_bc_to_device() { m_info->copy_bc_to_device(); }
215
217 [[nodiscard]] bool bc_initialized() const
218 {
219 return m_info->m_bc_copied_to_device;
220 }
221
223 [[nodiscard]] bool has_fillpatch_op() const
224 {
225 return static_cast<bool>(m_info->m_fillpatch_op);
226 }
227
236 const SimTime& time,
237 amrex::BCType::mathematicalBndryTypes bctype = amrex::BCType::hoextrap);
238
246 void to_uniform_space();
247
255 void to_stretched_space();
256
258 Field& state(FieldState fstate);
259 [[nodiscard]] const Field& state(FieldState fstate) const;
260
262 amrex::MultiFab& operator()(int lev);
263 const amrex::MultiFab& operator()(int lev) const;
264
266 amrex::Vector<amrex::MultiFab*> vec_ptrs();
267
269 [[nodiscard]] amrex::Vector<const amrex::MultiFab*> vec_const_ptrs() const;
270
272 void advance_states();
273
275 void copy_state(FieldState to_state, FieldState from_state);
276
279
280 // NOTE: We break naming rules for methods here to provide 1-1 mapping with
281 // the underlying MultiFab method names.
282
284 void setVal(amrex::Real value);
285
294 void
295 setVal(amrex::Real value, int start_comp, int num_comp = 1, int nghost = 0);
296
299 void setVal(const amrex::Vector<amrex::Real>& values, int nghost = 0);
300
303 ViewField<Field> subview(const int scomp = 0, const int ncomp = 1)
304 {
305 return {*this, scomp, ncomp};
306 }
307
312 template <typename T, class... Args>
313 void register_fill_patch_op(Args&&... args)
314 {
315 m_info->m_fillpatch_op.reset(new T(*this, std::forward<Args>(args)...));
316 }
317
320 template <typename T, class... Args>
321 void register_custom_bc(Args&&... args)
322 {
323 m_info->m_bc_func.emplace_back(
324 new T(*this, std::forward<Args>(args)...));
325 }
326
327 void fillpatch(amrex::Real time);
328 void fillpatch(amrex::Real time, amrex::IntVect ng);
330 amrex::Real time,
331 amrex::IntVect ng,
332 amrex::Array<Field*, AMREX_SPACEDIM>& fields) const;
333
334 void fillphysbc(amrex::Real time);
335 void fillphysbc(amrex::Real time, amrex::IntVect ng);
336
337 void apply_bc_funcs(FieldState rho_state);
338
339 void fillpatch(
340 int lev,
341 amrex::Real time,
342 amrex::MultiFab& mfab,
343 const amrex::IntVect& nghost);
344
346 int lev,
347 amrex::Real time,
348 amrex::MultiFab& mfab,
349 const amrex::IntVect& nghost);
350
351 void fillphysbc(
352 int lev,
353 amrex::Real time,
354 amrex::MultiFab& mfab,
355 const amrex::IntVect& nghost);
356
357 void set_inflow(
358 int lev,
359 amrex::Real time,
360 amrex::MultiFab& mfab,
361 const amrex::IntVect& nghost);
362
364 int lev,
365 amrex::Real time,
366 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM> mfabs);
367
369 const int lev,
370 const amrex::Real time,
371 amrex::MultiFab& mfab,
372 const int nghost)
373 {
374 fillpatch(lev, time, mfab, amrex::IntVect(nghost));
375 }
376
378 const int lev,
379 const amrex::Real time,
380 amrex::MultiFab& mfab,
381 const int nghost)
382 {
383 fillpatch_from_coarse(lev, time, mfab, amrex::IntVect(nghost));
384 }
385
387 const int lev,
388 const amrex::Real time,
389 amrex::MultiFab& mfab,
390 const int nghost)
391 {
392 fillphysbc(lev, time, mfab, amrex::IntVect(nghost));
393 }
394
396 const int lev,
397 const amrex::Real time,
398 amrex::MultiFab& mfab,
399 const int nghost)
400 {
401 set_inflow(lev, time, mfab, amrex::IntVect(nghost));
402 }
403
404 bool& in_uniform_space() { return m_mesh_mapped; }
405 [[nodiscard]] bool in_uniform_space() const { return m_mesh_mapped; }
406
408 [[nodiscard]] bool has_inout_bndry() const { return m_inout_bndry; }
409
412
413protected:
414 Field(
416 std::string name,
417 std::shared_ptr<FieldInfo> finfo,
418 unsigned fid,
420
423
425 const std::string m_name;
426
428 std::shared_ptr<FieldInfo> m_info;
429
431 const unsigned m_id;
432
435
439
441 bool m_mesh_mapped{false};
442
444 bool m_inout_bndry{false};
445};
446
447} // namespace amr_wind
448
449#endif /* FIELD_H */
Definition FieldBCOps.H:32
Definition FieldFillPatchOps.H:41
Definition Field.H:112
void set_default_fillpatch_bc(const SimTime &time, amrex::BCType::mathematicalBndryTypes bctype=amrex::BCType::hoextrap)
Definition Field.cpp:373
bool query_state(const FieldState fstate) const
Return true if the requested state exists for this field.
Definition Field.H:170
void fillphysbc(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const int nghost)
Definition Field.H:386
void to_stretched_space()
Definition Field.cpp:417
void to_uniform_space()
Definition Field.cpp:387
bool m_fillpatch_on_regrid
Flag indicating whether fill patch operation must be performed for this field during regrid.
Definition Field.H:438
void set_inout_bndry()
Set the inout_bndry flag.
Definition Field.H:411
void set_inflow(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const int nghost)
Definition Field.H:395
void fillpatch(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const int nghost)
Definition Field.H:368
bool has_inout_bndry() const
Check if any of the boundaries is a mass-inflow-outflow.
Definition Field.H:408
const amrex::Gpu::DeviceVector< amrex::BCRec > & bcrec_device() const
Definition Field.H:205
Field & state(FieldState fstate)
Return field at a different time state.
Definition Field.cpp:117
ViewField< Field > subview(const int scomp=0, const int ncomp=1)
Definition Field.H:303
amrex::Vector< amrex::MultiFab * > vec_ptrs()
Return a vector of MultiFab pointers for all levels.
Definition Field.cpp:145
const amrex::IntVect & num_grow() const
Ghost cells.
Definition Field.H:136
int num_time_states() const
Number of exact time states available for this field.
Definition Field.H:142
FieldLoc field_location() const
Cell, node, or face centered field.
Definition Field.H:155
unsigned id() const
Unique integer identifier for this field.
Definition Field.H:130
const std::string & base_name() const
Base name (without state information)
Definition Field.H:124
void set_inflow_sibling_fields(int lev, amrex::Real time, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > mfabs)
Definition Field.cpp:288
const amrex::GpuArray< const amrex::Real *, AMREX_SPACEDIM *2 > & bc_values_device() const
Definition Field.H:199
amrex::Vector< const amrex::MultiFab * > vec_const_ptrs() const
Return vector of const MultiFab* for all levels.
Definition Field.cpp:156
std::shared_ptr< FieldInfo > m_info
Common data for all field states.
Definition Field.H:428
bool & fillpatch_on_regrid()
Definition Field.H:163
FieldRepo & repo() const
FieldRepo instance that manages this field.
Definition Field.H:161
void advance_states()
Advance timestep for fields with multiple states.
Definition Field.cpp:300
bool & in_uniform_space()
Definition Field.H:404
void copy_state(FieldState to_state, FieldState from_state)
Copy a user-specified "from_state" to "to_state".
Definition Field.cpp:319
const amrex::GpuArray< BC, AMREX_SPACEDIM *2 > & bc_type() const
Definition Field.H:176
int num_states() const
Number of states available for this field.
Definition Field.H:145
Field(const Field &)=delete
void register_fill_patch_op(Args &&... args)
Definition Field.H:313
friend class FieldRepo
Definition Field.H:114
void fillphysbc(amrex::Real time)
Definition Field.cpp:265
const FieldState m_state
State for this field.
Definition Field.H:434
FieldState field_state() const
State of this field instance.
Definition Field.H:158
bool fillpatch_on_regrid() const
Definition Field.H:164
amrex::GpuArray< BC, AMREX_SPACEDIM *2 > & bc_type()
Definition Field.H:181
void register_custom_bc(Args &&... args)
Definition Field.H:321
void fillpatch_from_coarse(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const int nghost)
Definition Field.H:377
void copy_bc_to_device()
Copy BC data from host to device.
Definition Field.H:214
amrex::Vector< amrex::Vector< amrex::Real > > & bc_values()
Return reference to the host view of BC values array.
Definition Field.H:187
void fillpatch_from_coarse(int lev, amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost)
Definition Field.cpp:183
const std::string & name() const
Name of this field (including state information)
Definition Field.H:121
const unsigned m_id
Unique ID for this field.
Definition Field.H:431
bool bc_initialized() const
Return flag indicating whether BCs have been initialized.
Definition Field.H:217
void fillpatch_sibling_fields(amrex::Real time, amrex::IntVect ng, amrex::Array< Field *, AMREX_SPACEDIM > &fields) const
Definition Field.cpp:212
void apply_bc_funcs(FieldState rho_state)
Definition Field.cpp:267
int num_comp() const
Number of components for this field.
Definition Field.H:133
bool in_uniform_space() const
Definition Field.H:405
bool m_inout_bndry
Flag to indicate whether any of the boundaries is mass-inflow-outflow.
Definition Field.H:444
bool m_mesh_mapped
Flag to track mesh mapping (to uniform space) of field.
Definition Field.H:441
FieldRepo & m_repo
Reference to the FieldRepository instance.
Definition Field.H:422
const std::string m_name
Name of the field.
Definition Field.H:425
amrex::MultiFab & operator()(int lev)
Return MultiFab instance for a given level.
Definition Field.cpp:133
amrex::Vector< amrex::BCRec > & bcrec() const
Return reference to host view of BCRec array.
Definition Field.H:193
void set_inflow(int lev, amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost)
Definition Field.cpp:275
Field & operator=(const Field &)=delete
void setVal(amrex::Real value)
Set the field to a constant value at all levels.
Definition Field.cpp:341
bool has_fillpatch_op() const
Return a flag indicating whether a fillpatch Op has been registered.
Definition Field.H:223
Field & create_state(FieldState fstate)
Create a new time state after the field has been created.
Definition Field.cpp:331
void fillpatch(amrex::Real time)
Definition Field.cpp:210
Definition FieldRepo.H:86
Definition SimTime.H:33
Definition ViewField.H:24
FieldState
Definition FieldDescTypes.H:16
FieldLoc
Definition FieldDescTypes.H:29
This test case is intended as an evaluation of the momentum advection scheme.
Definition BCInterface.cpp:10
FieldLoc m_floc
Cell, node, face centered field type.
Definition Field.H:57
FieldInfo(std::string basename, int ncomp, int ngrow, int nstates, FieldLoc floc)
Definition Field.cpp:15
amrex::IntVect m_ngrow
Ghost cells.
Definition Field.H:51
void copy_bc_to_device()
Copy the BC information to device data structures.
Definition Field.cpp:65
amrex::Vector< amrex::Vector< amrex::Real > > m_bc_values
Definition Field.H:62
std::string m_basename
Field name without state information.
Definition Field.H:45
int m_nstates
Number of states available for this field.
Definition Field.H:54
std::unique_ptr< FieldFillPatchOpsBase > m_fillpatch_op
Function that handles filling patch and physics BC data for this field.
Definition Field.H:74
bool m_bc_copied_to_device
Flag indicating whether BCs have been initialized and copied to device.
Definition Field.H:80
amrex::GpuArray< BC, AMREX_SPACEDIM *2 > m_bc_type
Definition Field.H:61
amrex::GpuArray< const amrex::Real *, AMREX_SPACEDIM *2 > m_bc_values_d
Definition Field.H:63
bool bc_initialized()
Check indicating whether BCs are being set to sane values.
Definition Field.cpp:42
amrex::Vector< Field * > m_states
Vector holding references to available states for this field.
Definition Field.H:71
static constexpr int max_field_states
Maximum number of states allowed for a field.
Definition Field.H:31
amrex::Vector< amrex::BCRec > m_bcrec
Definition Field.H:66
amrex::Vector< std::unique_ptr< FieldBCIface > > m_bc_func
Custom boundary condition actions for this field.
Definition Field.H:77
amrex::Gpu::DeviceVector< amrex::Real > m_bc_values_dview
Definition Field.H:64
int m_ncomp
Number of components (Scalar = 1, Vector = 2, etc.)
Definition Field.H:48
amrex::Gpu::DeviceVector< amrex::BCRec > m_bcrec_d
Definition Field.H:67