/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
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 const int ncomp = 1,
142 const int ngrow = 0,
143 const int nstates = 1,
144 const 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 inline 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,
243 const FieldState fstate = FieldState::New) const;
244
247 Field& get_field(const int field_id) const
248 {
249 return *m_field_vec.at(field_id);
250 }
251
255
259
262 bool field_exists(
263 const std::string& name,
264 const FieldState fstate = FieldState::New) const;
265
277 const std::string& name,
278 const int ncomp = 1,
279 const int ngrow = 0,
280 const int nstates = 1,
281 const FieldLoc floc = FieldLoc::CELL);
282
285 const std::string& name,
286 const FieldState fstate = FieldState::New) const;
287
291 IntField& get_int_field(const int field_id) const
292 {
293 return *m_int_field_vec.at(field_id);
294 }
295
297 bool int_field_exists(
298 const std::string& name,
299 const FieldState fstate = FieldState::New) const;
300
311 std::unique_ptr<ScratchField> create_scratch_field(
312 const std::string& name,
313 const int ncomp = 1,
314 const int nghost = 0,
315 const FieldLoc floc = FieldLoc::CELL) const;
316
327 std::unique_ptr<ScratchField> create_scratch_field(
328 const int ncomp = 1,
329 const int nghost = 0,
330 const FieldLoc floc = FieldLoc::CELL) const;
331
332 std::unique_ptr<ScratchField> create_scratch_field_on_host(
333 const std::string& name,
334 const int ncomp = 1,
335 const int nghost = 0,
336 const FieldLoc floc = FieldLoc::CELL) const;
337
338 std::unique_ptr<ScratchField> create_scratch_field_on_host(
339 const int ncomp = 1,
340 const int nghost = 0,
341 const FieldLoc floc = FieldLoc::CELL) const;
342
343 std::unique_ptr<IntScratchField> create_int_scratch_field_on_host(
344 const std::string& name,
345 const int ncomp = 1,
346 const int nghost = 0,
347 const FieldLoc floc = FieldLoc::CELL) const;
348
349 std::unique_ptr<IntScratchField> create_int_scratch_field_on_host(
350 const int ncomp = 1,
351 const int nghost = 0,
352 const FieldLoc floc = FieldLoc::CELL) const;
353
355 void advance_states() noexcept;
356
358 const amrex::AmrCore& mesh() const { return m_mesh; }
359
361 int num_active_levels() const noexcept { return m_mesh.finestLevel() + 1; }
362
364 int num_fields() const noexcept
365 {
366 return static_cast<int>(m_field_vec.size());
367 }
368
370 const amrex::Vector<std::unique_ptr<Field>>& fields() const
371 {
372 return m_field_vec;
373 }
374
376 inline const amrex::FabFactory<amrex::FArrayBox>&
377 factory(int lev) const noexcept
378 {
379 return *m_leveldata[lev]->m_factory;
380 }
381
382protected:
388 inline amrex::MultiFab&
389 get_multifab(const unsigned fid, const int lev) noexcept
390 {
391 BL_ASSERT(lev <= m_mesh.finestLevel());
392 return m_leveldata[lev]->m_mfabs[fid];
393 }
394
400 inline amrex::iMultiFab&
401 get_int_fab(const unsigned fid, const int lev) noexcept
402 {
403 BL_ASSERT(lev <= m_mesh.finestLevel());
404 return m_leveldata[lev]->m_int_fabs[fid];
405 }
406
408 Field& create_state(Field& field, const FieldState fstate);
409
412 int lev,
413 const Field& field,
414 LevelDataHolder& level_data,
415 const amrex::FabFactory<amrex::FArrayBox>& factory);
416
418 void allocate_field_data(Field& field);
419
422 const amrex::BoxArray& ba,
423 const amrex::DistributionMapping& dm,
424 LevelDataHolder& level_data,
425 const amrex::FabFactory<amrex::FArrayBox>& factory);
426
428 int lev, const IntField& field, LevelDataHolder& level_data);
429
430 void allocate_field_data(const IntField& field);
431
433 const amrex::BoxArray& ba,
434 const amrex::DistributionMapping& dm,
435 LevelDataHolder& level_data,
436 const amrex::FabFactory<amrex::IArrayBox>& factory);
437
439 const amrex::AmrCore& m_mesh;
440
443 amrex::Vector<std::unique_ptr<LevelDataHolder>> m_leveldata;
444
446 mutable amrex::Vector<std::unique_ptr<Field>> m_field_vec;
447
449 mutable amrex::Vector<std::unique_ptr<IntField>> m_int_field_vec;
450
452 std::unordered_map<std::string, size_t> m_fid_map;
453
455 std::unordered_map<std::string, size_t> m_int_fid_map;
456
458 bool m_is_initialized{false};
459};
460
461} // namespace amr_wind
462
463#endif /* FIELDREPO_H */
Definition Field.H:116
Definition FieldRepo.H:86
void clear_level(int lev)
Remove a level during regrid.
Definition FieldRepo.cpp:77
int num_fields() const noexcept
Number of fields registered in the database.
Definition FieldRepo.H:364
std::unique_ptr< ScratchField > create_scratch_field(const std::string &name, const int ncomp=1, const int nghost=0, const FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:302
IntField & get_int_field(const int field_id) const
Definition FieldRepo.H:291
Field & declare_xf_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:182
std::unique_ptr< IntScratchField > create_int_scratch_field_on_host(const std::string &name, const int ncomp=1, const int nghost=0, const FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:368
const amrex::AmrCore & m_mesh
Reference to the mesh instance.
Definition FieldRepo.H:439
Field & get_field(const int field_id) const
Definition FieldRepo.H:247
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:455
bool m_is_initialized
Flag indicating if mesh is available to allocate field data.
Definition FieldRepo.H:458
Field & declare_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1, const FieldLoc floc=FieldLoc::CELL)
Definition FieldRepo.cpp:83
amrex::Vector< std::unique_ptr< IntField > > m_int_field_vec
Reference to integer field instances identified by unique integer.
Definition FieldRepo.H:449
Field & get_field(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:149
FieldRepo & operator=(const FieldRepo &)=delete
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
amrex::MultiFab & get_multifab(const unsigned fid, const int lev) noexcept
Definition FieldRepo.H:389
amrex::Vector< std::unique_ptr< Field > > m_field_vec
References to field instances identified by unique integer.
Definition FieldRepo.H:446
FieldRepo(const amrex::AmrCore &mesh)
Definition FieldRepo.H:91
amrex::Vector< std::unique_ptr< LevelDataHolder > > m_leveldata
Definition FieldRepo.H:443
void make_new_level_from_scratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition FieldRepo.cpp:12
bool int_field_exists(const std::string &name, const FieldState fstate=FieldState::New) const
Query if an integer field exists.
Definition FieldRepo.cpp:293
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
const amrex::FabFactory< amrex::FArrayBox > & factory(int lev) const noexcept
Return factory instance at a given level.
Definition FieldRepo.H:377
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:190
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:452
Field & declare_yf_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:197
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
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:431
bool field_exists(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:215
IntField & declare_int_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1, const FieldLoc floc=FieldLoc::CELL)
Definition FieldRepo.cpp:223
amrex::iMultiFab & get_int_fab(const unsigned fid, const int lev) noexcept
Definition FieldRepo.H:401
Field & get_mesh_mapping_field(FieldLoc floc) const
Definition FieldRepo.cpp:165
void make_new_level_from_coarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition FieldRepo.cpp:29
const amrex::Vector< std::unique_ptr< Field > > & fields() const
Return list of fields registered (for unit-test purposes only)
Definition FieldRepo.H:370
Field & create_state(Field &field, const FieldState fstate)
Create a new state for a field.
Definition FieldRepo.cpp:496
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:358
void remake_level(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition FieldRepo.cpp:53
std::unique_ptr< ScratchField > create_scratch_field_on_host(const std::string &name, const int ncomp=1, const int nghost=0, const FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:334
void advance_states() noexcept
Advance all fields with more than one timestate to the new timestep.
Definition FieldRepo.cpp:401
Definition IntField.H:20
FieldLoc
Definition FieldDescTypes.H:27
FieldState
Definition FieldDescTypes.H:14
@ NODE
Node-centered (e.g., for pressure)
@ ZFACE
Face-centered in z-direction.
@ XFACE
Face-centered in x-direction (e.g., face normal velocity)
@ CELL
Cell-centered (default)
@ YFACE
Face-centered in y-direction.
@ New
Same as FieldState::NP1.
Definition BCInterface.cpp:7
Definition console_io.cpp:25
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:7
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