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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/core/FieldRepo.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
FieldRepo.H
Go to the documentation of this file.
1#ifndef FIELDREPO_H
2#define FIELDREPO_H
3
4#include <string>
5#include <unordered_map>
6
13
14#include "AMReX_AmrCore.H"
15#include "AMReX_MultiFab.H"
16#include "AMReX_iMultiFab.H"
17
18namespace amr_wind {
19
32
36{
38
40 amrex::Vector<amrex::MultiFab> m_mfabs;
42 std::unique_ptr<amrex::FabFactory<amrex::FArrayBox>> m_factory;
43
45 amrex::Vector<amrex::iMultiFab> m_int_fabs;
46 std::unique_ptr<amrex::FabFactory<amrex::IArrayBox>> m_int_fact;
47};
48
86{
87public:
88 friend class Field;
89 friend class IntField;
90
91 explicit FieldRepo(const amrex::AmrCore& mesh)
92 : m_mesh(mesh), m_leveldata(mesh.maxLevel() + 1)
93 {}
94
95 FieldRepo(const FieldRepo&) = delete;
96 FieldRepo& operator=(const FieldRepo&) = delete;
97 ~FieldRepo() = default;
98
102 int lev,
103 amrex::Real time,
104 const amrex::BoxArray& ba,
105 const amrex::DistributionMapping& dm);
106
113 int lev,
114 amrex::Real time,
115 const amrex::BoxArray& ba,
116 const amrex::DistributionMapping& dm);
117
122 void remake_level(
123 int lev,
124 amrex::Real time,
125 const amrex::BoxArray& ba,
126 const amrex::DistributionMapping& dm);
127
129 void clear_level(int lev);
130
140 const std::string& name,
141 int ncomp = 1,
142 int ngrow = 0,
143 int nstates = 1,
144 FieldLoc floc = FieldLoc::CELL);
145
153 const std::string& name,
154 const int ncomp = 1,
155 const int ngrow = 0,
156 const int nstates = 1)
157 {
158 return declare_field(name, ncomp, ngrow, nstates, FieldLoc::CELL);
159 }
160
168 const std::string& name,
169 const int ncomp = 1,
170 const int ngrow = 0,
171 const int nstates = 1)
172 {
173 return declare_field(name, ncomp, ngrow, nstates, FieldLoc::NODE);
174 }
175
183 const std::string& name,
184 const int ncomp = 1,
185 const int ngrow = 0,
186 const int nstates = 1)
187 {
188 return declare_field(name, ncomp, ngrow, nstates, FieldLoc::XFACE);
189 }
190
198 const std::string& name,
199 const int ncomp = 1,
200 const int ngrow = 0,
201 const int nstates = 1)
202 {
203 return declare_field(name, ncomp, ngrow, nstates, FieldLoc::YFACE);
204 }
205
213 const std::string& name,
214 const int ncomp = 1,
215 const int ngrow = 0,
216 const int nstates = 1)
217 {
218 return declare_field(name, ncomp, ngrow, nstates, FieldLoc::ZFACE);
219 }
220
225 amrex::Vector<Field*> declare_face_normal_field(
226 const amrex::Vector<std::string>& names,
227 const int ncomp = 1,
228 const int ngrow = 0,
229 const int nstates = 1)
230 {
231 AMREX_ASSERT(names.size() == AMREX_SPACEDIM);
232 return {
233 &declare_field(names[0], ncomp, ngrow, nstates, FieldLoc::XFACE),
234 &declare_field(names[1], ncomp, ngrow, nstates, FieldLoc::YFACE),
235 &declare_field(names[2], ncomp, ngrow, nstates, FieldLoc::ZFACE),
236 };
237 }
238
242 const std::string& name, FieldState fstate = FieldState::New) const;
243
246 Field& get_field(const int field_id) const
247 {
248 return *m_field_vec.at(field_id);
249 }
250
254
258
261 bool field_exists(
262 const std::string& name, FieldState fstate = FieldState::New) const;
263
275 const std::string& name,
276 int ncomp = 1,
277 int ngrow = 0,
278 int nstates = 1,
279 FieldLoc floc = FieldLoc::CELL);
280
283 const std::string& name, FieldState fstate = FieldState::New) const;
284
288 IntField& get_int_field(const int field_id) const
289 {
290 return *m_int_field_vec.at(field_id);
291 }
292
294 bool int_field_exists(
295 const std::string& name, FieldState fstate = FieldState::New) const;
296
307 std::unique_ptr<ScratchField> create_scratch_field(
308 const std::string& name,
309 int ncomp = 1,
310 int nghost = 0,
311 FieldLoc floc = FieldLoc::CELL) const;
312
323 std::unique_ptr<ScratchField> create_scratch_field(
324 int ncomp = 1, int nghost = 0, FieldLoc floc = FieldLoc::CELL) const;
325
326 std::unique_ptr<ScratchField> create_scratch_field_on_host(
327 const std::string& name,
328 int ncomp = 1,
329 int nghost = 0,
330 FieldLoc floc = FieldLoc::CELL) const;
331
332 std::unique_ptr<ScratchField> create_scratch_field_on_host(
333 int ncomp = 1, int nghost = 0, FieldLoc floc = FieldLoc::CELL) const;
334
335 std::unique_ptr<IntScratchField> create_int_scratch_field_on_host(
336 const std::string& name,
337 int ncomp = 1,
338 int nghost = 0,
339 FieldLoc floc = FieldLoc::CELL) const;
340
341 std::unique_ptr<IntScratchField> create_int_scratch_field_on_host(
342 int ncomp = 1, int nghost = 0, FieldLoc floc = FieldLoc::CELL) const;
343
345 void advance_states();
346
348 const amrex::AmrCore& mesh() const { return m_mesh; }
349
351 int num_active_levels() const { return m_mesh.finestLevel() + 1; }
352
354 int num_fields() const { return static_cast<int>(m_field_vec.size()); }
355
357 const amrex::Vector<std::unique_ptr<Field>>& fields() const
358 {
359 return m_field_vec;
360 }
361
363 const amrex::FabFactory<amrex::FArrayBox>& factory(int lev) const
364 {
365 return *m_leveldata[lev]->m_factory;
366 }
367
368protected:
374 amrex::MultiFab& get_multifab(const unsigned fid, const int lev)
375 {
376 BL_ASSERT(lev <= m_mesh.finestLevel());
377 return m_leveldata[lev]->m_mfabs[fid];
378 }
379
385 amrex::iMultiFab& get_int_fab(const unsigned fid, const int lev)
386 {
387 BL_ASSERT(lev <= m_mesh.finestLevel());
388 return m_leveldata[lev]->m_int_fabs[fid];
389 }
390
392 Field& create_state(Field& field, FieldState fstate);
393
396 int lev,
397 const Field& field,
398 LevelDataHolder& level_data,
399 const amrex::FabFactory<amrex::FArrayBox>& factory);
400
402 void allocate_field_data(Field& field);
403
406 const amrex::BoxArray& ba,
407 const amrex::DistributionMapping& dm,
408 LevelDataHolder& level_data,
409 const amrex::FabFactory<amrex::FArrayBox>& factory);
410
412 int lev, const IntField& field, LevelDataHolder& level_data);
413
414 void allocate_field_data(const IntField& field);
415
417 const amrex::BoxArray& ba,
418 const amrex::DistributionMapping& dm,
419 LevelDataHolder& level_data,
420 const amrex::FabFactory<amrex::IArrayBox>& factory);
421
423 const amrex::AmrCore& m_mesh;
424
427 amrex::Vector<std::unique_ptr<LevelDataHolder>> m_leveldata;
428
430 mutable amrex::Vector<std::unique_ptr<Field>> m_field_vec;
431
433 mutable amrex::Vector<std::unique_ptr<IntField>> m_int_field_vec;
434
436 std::unordered_map<std::string, size_t> m_fid_map;
437
439 std::unordered_map<std::string, size_t> m_int_fid_map;
440
442 bool m_is_initialized{false};
443};
444
445} // namespace amr_wind
446
447#endif /* FIELDREPO_H */
Definition Field.H:112
void clear_level(int lev)
Remove a level during regrid.
Definition FieldRepo.cpp:80
IntField & get_int_field(const int field_id) const
Definition FieldRepo.H:288
std::unique_ptr< ScratchField > create_scratch_field_on_host(const std::string &name, int ncomp=1, int nghost=0, FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:337
Field & declare_xf_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:182
amrex::MultiFab & get_multifab(const unsigned fid, const int lev)
Definition FieldRepo.H:374
const amrex::AmrCore & m_mesh
Reference to the mesh instance.
Definition FieldRepo.H:423
Field & get_field(const int field_id) const
Definition FieldRepo.H:246
std::unordered_map< std::string, size_t > m_int_fid_map
Map of integer field name to unique integer ID for lookups.
Definition FieldRepo.H:439
bool m_is_initialized
Flag indicating if mesh is available to allocate field data.
Definition FieldRepo.H:442
IntField & declare_int_field(const std::string &name, int ncomp=1, int ngrow=0, int nstates=1, FieldLoc floc=FieldLoc::CELL)
Definition FieldRepo.cpp:226
const amrex::FabFactory< amrex::FArrayBox > & factory(int lev) const
Return factory instance at a given level.
Definition FieldRepo.H:363
amrex::Vector< std::unique_ptr< IntField > > m_int_field_vec
Reference to integer field instances identified by unique integer.
Definition FieldRepo.H:433
amrex::iMultiFab & get_int_fab(const unsigned fid, const int lev)
Definition FieldRepo.H:385
Field & create_state(Field &field, FieldState fstate)
Create a new state for a field.
Definition FieldRepo.cpp:499
bool int_field_exists(const std::string &name, FieldState fstate=FieldState::New) const
Query if an integer field exists.
Definition FieldRepo.cpp:296
FieldRepo & operator=(const FieldRepo &)=delete
amrex::Vector< std::unique_ptr< Field > > m_field_vec
References to field instances identified by unique integer.
Definition FieldRepo.H:430
FieldRepo(const amrex::AmrCore &mesh)
Definition FieldRepo.H:91
IntField & get_int_field(const std::string &name, FieldState fstate=FieldState::New) const
Return a reference to an integer field.
Definition FieldRepo.cpp:279
void advance_states()
Advance all fields with more than one timestate to the new timestep.
Definition FieldRepo.cpp:404
amrex::Vector< std::unique_ptr< LevelDataHolder > > m_leveldata
Array (size: nlevels) of data holder that contains another array of MultiFabs for all fields at that ...
Definition FieldRepo.H:427
void make_new_level_from_scratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Perform field data management tasks when amrex::AmrCore::MakeNewLevelFromScratch is called.
Definition FieldRepo.cpp:15
Field & declare_nd_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:167
amrex::Vector< Field * > declare_face_normal_field(const amrex::Vector< std::string > &names, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:225
Field & declare_zf_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:212
Field & get_mesh_mapping_det_j(FieldLoc floc) const
Definition FieldRepo.cpp:193
FieldRepo(const FieldRepo &)=delete
std::unordered_map< std::string, size_t > m_fid_map
Map of field name to unique integer ID for lookups.
Definition FieldRepo.H:436
Field & declare_yf_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:197
std::unique_ptr< IntScratchField > create_int_scratch_field_on_host(const std::string &name, int ncomp=1, int nghost=0, FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:371
bool field_exists(const std::string &name, FieldState fstate=FieldState::New) const
Query if field uniquely identified by name and time state exists in repository.
Definition FieldRepo.cpp:218
void allocate_field_data(int lev, const Field &field, LevelDataHolder &level_data, const amrex::FabFactory< amrex::FArrayBox > &factory)
Allocate field data for a single level outside of regrid.
Definition FieldRepo.cpp:434
int num_fields() const
Number of fields registered in the database.
Definition FieldRepo.H:354
friend class Field
Definition FieldRepo.H:88
Field & declare_field(const std::string &name, int ncomp=1, int ngrow=0, int nstates=1, FieldLoc floc=FieldLoc::CELL)
Definition FieldRepo.cpp:86
int num_active_levels() const
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:351
Field & get_mesh_mapping_field(FieldLoc floc) const
Definition FieldRepo.cpp:168
void make_new_level_from_coarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition FieldRepo.cpp:32
Field & get_field(const std::string &name, FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:152
const amrex::Vector< std::unique_ptr< Field > > & fields() const
Return list of fields registered (for unit-test purposes only)
Definition FieldRepo.H:357
std::unique_ptr< ScratchField > create_scratch_field(const std::string &name, int ncomp=1, int nghost=0, FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:305
friend class IntField
Definition FieldRepo.H:89
Field & declare_cc_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:152
const amrex::AmrCore & mesh() const
Return a reference to the underlying AMR mesh instance.
Definition FieldRepo.H:348
void remake_level(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition FieldRepo.cpp:56
Definition IntField.H:20
FieldState
Definition FieldDescTypes.H:16
FieldLoc
Definition FieldDescTypes.H:29
@ New
Same as FieldState::NP1.
Definition FieldDescTypes.H:22
@ NODE
Node-centered (e.g., for pressure)
Definition FieldDescTypes.H:31
@ ZFACE
Face-centered in z-direction.
Definition FieldDescTypes.H:34
@ XFACE
Face-centered in x-direction (e.g., face normal velocity)
Definition FieldDescTypes.H:32
@ CELL
Cell-centered (default)
Definition FieldDescTypes.H:30
@ YFACE
Face-centered in y-direction.
Definition FieldDescTypes.H:33
This test case is intended as an evaluation of the momentum advection scheme.
Definition BCInterface.cpp:10
Definition FieldRepo.H:36
amrex::Vector< amrex::iMultiFab > m_int_fabs
int fabs for all known fields at this level
Definition FieldRepo.H:45
std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > m_factory
Factory for creating new FABs.
Definition FieldRepo.H:42
LevelDataHolder()
Definition FieldRepo.cpp:10
amrex::Vector< amrex::MultiFab > m_mfabs
real multifabs for all the known fields at this level
Definition FieldRepo.H:40
std::unique_ptr< amrex::FabFactory< amrex::IArrayBox > > m_int_fact
Definition FieldRepo.H:46