18#include <DisjointBoxLayout.H>
19#include <ProblemDomain.H>
35#include <CD_NamespaceHeader.H>
108 template <
typename T>
130 template <
typename T>
155 template <
typename T>
171 template <
typename T>
179 template <
typename T>
187 template <
typename T>
197 template <
typename T>
207 template <
typename T,
typename S>
216 template <
typename T>
225 template <
typename T>
234 template <
typename T>
244 template <
typename T>
254 template <
typename T>
265 template <
typename T>
373 const phase::which_phase
a_phase,
386 const phase::which_phase
a_phase,
457 const phase::which_phase
a_phase)
const;
472 const phase::which_phase
a_phase)
const;
508 const phase::which_phase
a_phase,
525 const phase::which_phase
a_phase,
542 const phase::which_phase
a_phase,
558 const phase::which_phase
a_phase,
573 const phase::which_phase
a_phase,
775 const phase::which_phase
a_phase,
789 const phase::which_phase
a_phase,
812 const phase::which_phase
a_phase,
834 const phase::which_phase
a_phase,
856 const phase::which_phase
a_phase,
869 const phase::which_phase
a_phase,
911 const phase::which_phase
a_phase,
953 template <
class P,
typename Ret, Ret (P::*MemberFunc)() const>
957 const phase::which_phase&
a_phase,
972 template <
class P,
typename Ret, Ret (P::*MemberFunc)() const>
976 const phase::which_phase&
a_phase,
990 template <
class P,
class Ret, Ret (P::*MemberFunc)()>
994 const phase::which_phase&
a_phase,
1012 const phase::which_phase&
a_phase,
1028 const phase::which_phase&
a_phase,
1058 const phase::which_phase&
a_phase,
1076 const phase::which_phase&
a_phase,
1093 const phase::which_phase&
a_phase)
const;
1107 const phase::which_phase
a_phase,
1137 const phase::which_phase
a_phase,
1165 const phase::which_phase
a_phase,
1194 const phase::which_phase
a_phase)
const;
1278 const phase::which_phase
a_phase,
1299 const phase::which_phase
a_phase,
1335 const phase::which_phase
a_phase,
1336 const int a_level)
const noexcept;
1351 const phase::which_phase
a_phase,
1367 const phase::which_phase
a_phase)
const noexcept;
1381 const phase::which_phase
a_phase,
1382 const int a_level)
const noexcept;
1395 const phase::which_phase
a_phase)
const noexcept;
1409 const phase::which_phase
a_phase,
1410 const int a_level)
const noexcept;
1425 const phase::which_phase
a_phase,
1442 const phase::which_phase&
a_phase)
const noexcept;
1459 const phase::which_phase&
a_phase)
const noexcept;
2194 const phase::which_phase
a_phase,
2195 const int a_lvl)
const;
2198#include <CD_NamespaceFooter.H>
Implementation of CD_AmrMesh.H.
Average
Various averaging methods.
Definition CD_Average.H:24
BoxSorting
Enum for sorting boxes.
Definition CD_BoxSorting.H:21
CoarseFineDeposition
Coarse-fine deposition types (see CD_EBAMRParticleMesh for how these are handled).
Definition CD_CoarseFineDeposition.H:26
Declaration of base class for defining geometries.
Simple enum for distinguishing copying strategies.
CopyStrategy
Enum for distinguishing how we copy data Valid => valid region ValidGhost => valid+ghost region.
Definition CD_CopyStrategy.H:23
DepositionType
Deposition types.
Definition CD_DepositionType.H:23
Declaration of factory class for DomainFluxIFFAB.
Declaration of a BaseIFFAB wrapper that holds domain fluxes.
Class for holding data across EBAMR hierarchies.
Declaration of conservative coarsening utility.
Declaration of a class which can reflux over the coarse-fine interface.
Declaration of a static class for various load balancing operations.
Declaration of a wrapper for wrapping multifluid EBLevelGrids.
Declaration of a class for holding particles on an AMR hierarchy.
Declaration of the Realm class.
Class for handling spatial operations.
Definition CD_AmrMesh.H:42
AmrMesh & operator=(const AmrMesh &a_other)=delete
Disallowed copy assignment.
void removeCoveredParticlesVoxels(ParticleContainer< P > &a_particles, const phase::which_phase &a_phase) const
Function which removes particles from the domain if they fall inside the EB.
Definition CD_AmrMeshImplem.H:605
void computeGradient(EBAMRCellData &a_gradient, const EBAMRCellData &a_phi, const std::string a_realm, const phase::which_phase a_phase) const
Compute cell-centered gradient over an AMR hierarchy.
Definition CD_AmrMesh.cpp:1334
void transferCoveredParticlesDiscrete(ParticleContainer< P > &a_particlesFrom, ParticleContainer< P > &a_particlesTo, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which transferse particles from one particle container to another if they fall inside the EB...
Definition CD_AmrMeshImplem.H:738
int m_refRatio
Refinement ratio.
Definition CD_AmrMesh.H:1936
void parseRuntimeOptions()
Parse runtime options.
Definition CD_AmrMesh.cpp:946
std::map< std::string, Vector< Copier > > m_oldToNewCellCopiers
Storage for copiers from the old grids to the new ones.
Definition CD_AmrMesh.H:1856
void allocatePointer(Vector< RefCountedPtr< T > > &a_data) const
Allocate pointer but not any memory blocks.
Definition CD_AmrMeshImplem.H:303
std::vector< std::string > getRealms() const
Get the name of all Realms.
Definition CD_AmrMesh.cpp:3661
RealVect getProbHi() const
Get upper-left corner of computational domain.
Definition CD_AmrMesh.cpp:2942
void parseBrBufferSize()
Parse buffer size for Berger-Rigoutsous grid algorithm.
Definition CD_AmrMesh.cpp:2635
Real getFinestDx() const
Get resolution on the finest grid level.
Definition CD_AmrMesh.cpp:3093
const AMRMask & getValidCells(const std::string a_realm) const
Get a map of all valid cells on a specified realm.
Definition CD_AmrMesh.cpp:3225
const RefCountedPtr< BaseIF > & getBaseImplicitFunction(const phase::which_phase a_phase) const
Get implicit function for a specific phase.
Definition CD_AmrMesh.cpp:3126
AmrMesh(const AmrMesh &&a_other)=delete
Disallowed move constructor.
void harmonicAverage(MFAMRCellData &a_data, const std::string a_realm) const
Harmonic coarsening of multifluid data.
Definition CD_AmrMesh.cpp:1520
Real m_fillRatioBR
Fill ratio.
Definition CD_AmrMesh.H:1921
ProblemDomain getFinestDomain() const
Get finest domain.
Definition CD_AmrMesh.cpp:3082
void parseCellCentroidInterpolation()
Parse centroid interpolation stencils.
Definition CD_AmrMesh.cpp:2786
int m_oldFinestLevel
Finest level before a regrid.
Definition CD_AmrMesh.H:1951
void copyData(EBAMRData< T > &a_dst, const EBAMRData< T > &a_src, const CopyStrategy &a_toRegion=CopyStrategy::Valid, const CopyStrategy &a_fromRegion=CopyStrategy::Valid) const noexcept
Method for copying from a source container to a destination container. User supplies information abou...
Definition CD_AmrMeshImplem.H:22
const Vector< DisjointBoxLayout > & getProxyGrids() const
Get the "proxy" grids in AmrMesh.
Definition CD_AmrMesh.cpp:3148
void transferCoveredParticlesIF(ParticleContainer< P > &a_particlesFrom, ParticleContainer< P > &a_particlesTo, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which transferse particles from one particle container to another if they fall inside the EB...
Definition CD_AmrMeshImplem.H:664
IntVect m_numCells
Coarsest box where we compute.
Definition CD_AmrMesh.H:1916
void postRegrid()
Run post-regrid operations.
Definition CD_AmrMesh.cpp:1039
void parseMaxAmrDepth()
Parse the maximum permitted AMR depth.
Definition CD_AmrMesh.cpp:2471
Vector< RefCountedPtr< EBMultigridInterpolator > > & getMultigridInterpolator(const std::string a_realm, const phase::which_phase a_phase) const
Get multigrid interpolation utility.
Definition CD_AmrMesh.cpp:3375
RealVect m_probHi
Domain simulation corner.
Definition CD_AmrMesh.H:1931
int m_maxBoxSize
Max box size.
Definition CD_AmrMesh.H:1966
void registerOperator(const std::string a_operator, const std::string a_realm, const phase::which_phase a_phase)
Register an operator over a realm and a phase.
Definition CD_AmrMesh.cpp:3541
Vector< DisjointBoxLayout > m_grids
Grids.
Definition CD_AmrMesh.H:2031
RefCountedPtr< MultiFluidIndexSpace > m_multifluidIndexSpace
MultiFluidIndexSpace.
Definition CD_AmrMesh.H:1901
RealVect m_probLo
Domain simulation corner.
Definition CD_AmrMesh.H:1926
void nonConservativeDivergence(EBAMRIVData &a_nonConsDivF, const EBAMRCellData &a_kappaDivF, const std::string &a_realm, const phase::which_phase &a_phase) const noexcept
Compute a non-conservative divergence.
Definition CD_AmrMesh.cpp:3461
void remapToNewGrids(ParticleContainer< P > &a_particles, const int a_lmin, const int a_newFinestLevel) const noexcept
Regrid particle to new grids.
Definition CD_AmrMeshImplem.H:362
int m_bufferSizeBR
Set buffer size.
Definition CD_AmrMesh.H:1976
void setCoarsestGrid(const IntVect &a_nCells)
Set the coarsest grid cells.
Definition CD_AmrMesh.cpp:2430
void parseMaxSimulationDepth()
Set maximum simulation depth.
Definition CD_AmrMesh.cpp:2490
int getMaxSimulationDepth() const
Get maximum permitted simulation depth.
Definition CD_AmrMesh.cpp:2975
void interpGhostPwl(MFAMRCellData &a_data, const std::string a_realm) const
Interpolate ghost cells over a realm and phase. This uses piecewise linear interpolation (with limite...
Definition CD_AmrMesh.cpp:1950
int m_maxSimulationDepth
Maximum allowed depth for simulation.
Definition CD_AmrMesh.H:1961
void setFinestLevel(const int a_finestLevel)
Set the finest level.
Definition CD_AmrMesh.cpp:2547
void transferIrregularParticles(ParticleContainer< P > &a_dstParticles, ParticleContainer< P > &a_srcParticles, const phase::which_phase a_phase, const std::function< void(P &)> a_transferModifier=[](P &) -> void { return;}) const noexcept
Transfer particles that are on the wrong side of the EB to a different container.
Definition CD_AmrMeshImplem.H:897
bool queryRealm(const std::string a_realm) const
Query if a realm exists.
Definition CD_AmrMesh.cpp:3500
virtual ~AmrMesh()
Destructor.
Definition CD_AmrMesh.cpp:46
void removeCoveredParticlesIF(ParticleContainer< P > &a_particles, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which removes particles from the domain if they fall inside the EB.
Definition CD_AmrMeshImplem.H:449
int m_multigridInterpRadius
Multigrid interpolation radius.
Definition CD_AmrMesh.H:2006
void arithmeticAverage(MFAMRCellData &a_data, const std::string a_realm) const
Arithmetic coarsening of multifluid data.
Definition CD_AmrMesh.cpp:1509
void buildDomains()
Build domains.
Definition CD_AmrMesh.cpp:982
GridGenerationMethod m_gridGenerationMethod
Grid generation method.
Definition CD_AmrMesh.H:1891
void interpGhostMG(MFAMRCellData &a_data, const std::string a_realm) const
Interpolate ghost cells over a realm and phase.
Definition CD_AmrMesh.cpp:2013
void parseRefinementRatios()
Parse refinement ratios.
Definition CD_AmrMesh.cpp:2509
Vector< int > m_refinementRatios
AMR resolutions.
Definition CD_AmrMesh.H:2041
int m_verbosity
Verbosity.
Definition CD_AmrMesh.H:1941
void conservativeAverage(MFAMRCellData &a_data, const std::string a_realm) const
Conservative coarsening of multifluid data.
Definition CD_AmrMesh.cpp:1531
bool m_hasRegridCopiers
Has regrid copiers or not.
Definition CD_AmrMesh.H:2026
Vector< RefCountedPtr< EBCoarseToFineInterp > > & getFineInterp(const std::string a_realm, const phase::which_phase a_phase) const
Get interpolator.
Definition CD_AmrMesh.cpp:3392
void parseMaxBoxSize()
Parse the maximum permitted box size.
Definition CD_AmrMesh.cpp:2597
int m_multigridInterpOrder
Multigrid interpolation order.
Definition CD_AmrMesh.H:2001
std::map< std::pair< std::string, std::string >, Vector< Copier > > m_validToValidRealmCopiers
Map for copying between various Realms. First index is the "from" realm and second index is the "to" ...
Definition CD_AmrMesh.H:1868
std::map< std::pair< std::string, std::string >, Vector< Copier > > m_validGhostToValidRealmCopiers
Map for copying between various Realms. First index is the "from" realm and second index is the "to" ...
Definition CD_AmrMesh.H:1880
void intersectParticlesRaycastIF(ParticleContainer< P > &a_activeParticles, ParticleContainer< P > &a_ebParticles, ParticleContainer< P > &a_domainParticles, const phase::which_phase a_phase, const Real a_tolerance, const bool a_deleteParticles, const std::function< void(P &)> a_nonDeletionModifier=[](P &) -> void { return;}) const noexcept
Particle intersection algorithm based on ray-casting.
Definition CD_AmrMeshImplem.H:981
void parseNumGhostCells()
Parse the number of ghost cells.
Definition CD_AmrMesh.cpp:2728
void sanityCheck() const
Do a sanity check to make sure everything is set up correctly.
Definition CD_AmrMesh.cpp:2872
int m_numGhostCells
Number of ghost cells.
Definition CD_AmrMesh.H:1991
void buildGrids(const Vector< IntVectSet > &a_tags, const int a_lmin, const int a_hardcap=-1)
Build new internal AMR grids.
Definition CD_AmrMesh.cpp:1128
void parseVerbosity()
Parse the verbosity for AmrMesh.
Definition CD_AmrMesh.cpp:2441
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGrid(const std::string a_realm, const phase::which_phase a_phase) const
Get the EBLevelGrid for a Realm and phase.
Definition CD_AmrMesh.cpp:3257
GridGenerationMethod
Enum for having understandable notation for grid generation.
Definition CD_AmrMesh.H:1832
std::map< std::pair< std::string, std::string >, Vector< Copier > > m_validGhostToValidGhostRealmCopiers
Map for copying between various Realms. First index is the "from" realm and second index is the "to" ...
Definition CD_AmrMesh.H:1886
void registerRealm(const std::string a_realm)
Register a new realm.
Definition CD_AmrMesh.cpp:3528
void removeCoveredParticlesDiscrete(ParticleContainer< P > &a_particles, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which removes particles from the domain if they fall inside the EB.
Definition CD_AmrMeshImplem.H:515
bool getEbCf() const
Check if mesh has an EBCF.
Definition CD_AmrMesh.cpp:3517
Vector< Real > m_dx
Level resolutions.
Definition CD_AmrMesh.H:2046
const Vector< DisjointBoxLayout > & getGrids(const std::string a_realm) const
Get the grids.
Definition CD_AmrMesh.cpp:3159
void intersectParticlesBisectIF(ParticleContainer< P > &a_activeParticles, ParticleContainer< P > &a_ebParticles, ParticleContainer< P > &a_domainParticles, const phase::which_phase a_phase, const Real a_bisectionStep, const bool a_deleteParticles, const std::function< void(P &)> a_nonDeletionModifier=[](P &) -> void { return;}) const noexcept
Particle intersection algorithm based on bisection.
Definition CD_AmrMeshImplem.H:1153
int m_numEbGhostsCells
Number of ghost cells to use for eb stuff.
Definition CD_AmrMesh.H:1986
void parseOptions()
Parse options. Called during the constructor.
Definition CD_AmrMesh.cpp:912
void buildCopiers()
Build copiers for copying between realms.
Definition CD_AmrMesh.cpp:1269
RealVect getProbLo() const
Get lower-left corner of computational domain.
Definition CD_AmrMesh.cpp:2931
BoxSorting m_boxSort
Box sorting.
Definition CD_AmrMesh.H:1896
void setMultifluidIndexSpace(const RefCountedPtr< MultiFluidIndexSpace > &a_multiFluidIndexSpace)
Sets multifluid index space.
Definition CD_AmrMesh.cpp:888
int m_maxAmrDepth
Maximum amr depth.
Definition CD_AmrMesh.H:1956
void defineRealms()
Define Realms.
Definition CD_AmrMesh.cpp:3574
int m_finestLevel
Finest level.
Definition CD_AmrMesh.H:1946
void reallocate(EBAMRCellData &a_data, const phase::which_phase a_phase, const int a_lmin) const
Reallocate data.
Definition CD_AmrMesh.cpp:582
std::map< std::string, Vector< Copier > > m_oldToNewEBCopiers
Storage for copiers from the old grids to the new ones.
Definition CD_AmrMesh.H:1862
const Vector< EBISLayout > & getEBISLayout(const std::string a_realm, const phase::which_phase a_phase) const
Get EBISLayouts for a Realm and phase.
Definition CD_AmrMesh.cpp:3175
Vector< RefCountedPtr< EBReflux > > & getFluxRegister(const std::string a_realm, const phase::which_phase a_phase) const
Get flux register.
Definition CD_AmrMesh.cpp:3409
void depositParticles(EBAMRCellData &a_meshData, const std::string &a_realm, const phase::which_phase &a_phase, const ParticleContainer< P > &a_particles, const DepositionType a_depositionType, const CoarseFineDeposition a_coarseFineDeposition, const bool a_forceIrregNGP=false)
Deposit scalar particle quantities on the mesh.
Definition CD_AmrMeshImplem.H:383
Vector< RefCountedPtr< EBFluxRedistribution > > & getRedistributionOp(const std::string a_realm, const phase::which_phase a_phase) const
Get the redistribution operators.
Definition CD_AmrMesh.cpp:3426
std::map< phase::which_phase, RefCountedPtr< BaseIF > > m_baseif
Implicit functions.
Definition CD_AmrMesh.H:1845
int getNumberOfGhostCells() const
Get the default number of ghost cells.
Definition CD_AmrMesh.cpp:2986
void interpGhost(EBAMRCellData &a_data, const std::string a_realm, const phase::which_phase a_phase) const
Interpolate ghost vectors over a realm, using the default ghost cell interpolation method.
Definition CD_AmrMesh.cpp:1860
void interpToCentroids(EBAMRCellData &a_data, const std::string a_realm, const phase::which_phase a_phase) const noexcept
Interpolate data to centroids on realm and phase.
Definition CD_AmrMesh.cpp:2252
const Vector< RefCountedPtr< MFLevelGrid > > & getMFLevelGrid(const std::string a_realm) const
Get EBISLayouts for a Realm.
Definition CD_AmrMesh.cpp:3291
void regridAmr(const Vector< IntVectSet > &a_tags, const int a_lmin, const int a_hardcap=-1)
Regrid AMR. This versions generates the grids and Realms, but not the operator.
Definition CD_AmrMesh.cpp:1076
AmrMesh & operator=(const AmrMesh &&a_other)=delete
Disallowed move assignment.
void parseCoarsestLevelNumCells()
Parse the coarsest domain grid.
Definition CD_AmrMesh.cpp:2454
RefCountedPtr< EBIndexSpace > & getEBIndexSpace(const phase::which_phase a_phase)
Get EBIndexSpace corresponding to a particular phase.
Definition CD_AmrMesh.H:1822
void parseRedistributionRadius()
Parse the default redistribution radius.
Definition CD_AmrMesh.cpp:2770
void allocate(Vector< RefCountedPtr< ParticleData< T > > > &a_particles, const std::string a_realm) const
Template class for generic allocation of particle data.
Definition CD_AmrMeshImplem.H:237
void setGrids(const Vector< Vector< Box > > &a_boxes, const std::map< std::string, Vector< Vector< long int > > > &a_realmsAndLoads)
Set grids from boxes and computational loads.
Definition CD_AmrMesh.cpp:2560
std::map< std::pair< std::string, std::string >, Vector< Copier > > m_validToValidGhostRealmCopiers
Map for copying between various Realms. First index is the "from" realm and second index is the "to" ...
Definition CD_AmrMesh.H:1874
const Vector< Real > & getDx() const
Get spatial resolutions.
Definition CD_AmrMesh.cpp:3104
int getBrBuffer() const
Return buffer for B-R mesh refinement algorithm.
Definition CD_AmrMesh.cpp:3041
int getMaxBoxSize() const
Get maximum permitted box size.
Definition CD_AmrMesh.cpp:3030
std::map< std::string, RefCountedPtr< Realm > > m_realms
These are all the Realms.
Definition CD_AmrMesh.H:1840
void interpToEB(EBAMRIVData &a_centroidData, const EBAMRCellData &a_cellData, const std::string a_realm, const phase::which_phase a_phase) const noexcept
Interpolate data to EB centroids on realm and phase.
Definition CD_AmrMesh.cpp:2361
int getMaxEbisBoxSize() const
Get maximum box size for EBIS generation.
Definition CD_AmrMesh.cpp:3052
AmrMesh()
Default constructor.
Definition CD_AmrMesh.cpp:31
void parseEbGhostCells()
Parse number of ghost cells for eb stuff.
Definition CD_AmrMesh.cpp:2712
void parseMultigridInterpolator()
Parse settings for the multigrid interpolator.
Definition CD_AmrMesh.cpp:2748
Vector< RefCountedPtr< EBCoarAve > > & getCoarseAverage(const std::string a_realm, const phase::which_phase a_phase) const
Get the coarsening utility.
Definition CD_AmrMesh.cpp:3358
EBAMRCellData slice(EBAMRCellData &a_original, const Interval a_variables) const noexcept
Slice cell-centered data in order to fetch a subset of components.
Definition CD_AmrMesh.cpp:50
int m_blockingFactor
Blocking factor.
Definition CD_AmrMesh.H:1981
int getRefinementRatio(const int a_level1, const int a_level2) const
Get refinement factor between two levels.
Definition CD_AmrMesh.cpp:3063
void interpolateParticles(ParticleContainer< P > &a_particles, const std::string &a_realm, const phase::which_phase &a_phase, const EBAMRCellData &a_meshScalarField, const DepositionType a_interpType, const bool a_forceIrregNGP=false) const
Interpolate mesh data onto a particle position.
Definition CD_AmrMeshImplem.H:427
void parseBlockingFactor()
Parse the B-R blocking factor. For tiled mesh refinement this parses the tile size.
Definition CD_AmrMesh.cpp:2696
int getBlockingFactor() const
Get blocking factor (i.e. the smallest possible box).
Definition CD_AmrMesh.cpp:3019
CellCentroidInterpolation::Type m_cellCentroidInterpolationType
Interpolation type for cell to cell centroid interpolation.
Definition CD_AmrMesh.H:1906
void parseEBCentroidInterpolation()
Parse EB interpolation (or extrapolation) stencils.
Definition CD_AmrMesh.cpp:2829
int m_multigridInterpWeight
Multigrid interpolation weights.
Definition CD_AmrMesh.H:2011
std::map< std::string, Vector< DisjointBoxLayout > > m_oldGrids
Old grids.
Definition CD_AmrMesh.H:1850
EBCentroidInterpolation::Type m_ebCentroidInterpolationType
Interpolation type for cell to EB centroid interpolation
Definition CD_AmrMesh.H:1911
const AMRMask & getMask(const std::string a_mask, const int a_buffer, const std::string a_realm) const
Get a registered mask.
Definition CD_AmrMesh.cpp:3209
void average(MFAMRCellData &a_data, const std::string a_realm, const Average &a_average) const
Average multifluid data over a specified realm.
Definition CD_AmrMesh.cpp:1469
EBAMRSurfaceDeposition & getSurfaceDeposition(const std::string a_realm, const phase::which_phase a_phase) const
Get EBAMRSurfaceDeposition surface deposition operator.
Definition CD_AmrMesh.cpp:3341
void preRegrid()
Run pre-regrid operations.
Definition CD_AmrMesh.cpp:1007
void registerMask(const std::string a_mask, const int a_buffer, const std::string a_realm)
Register a boolean mask over a realm.
Definition CD_AmrMesh.cpp:3558
AmrMesh(const AmrMesh &a_other)=delete
Disallowed copy constructor.
void parseBrFillRatio()
Parse the Berger-Rigoutsos fill ratio.
Definition CD_AmrMesh.cpp:2531
Vector< ProblemDomain > m_domains
Problem domains.
Definition CD_AmrMesh.H:2036
void parseProbLoHiCorners()
Parse the low/high corners of the computational domain.
Definition CD_AmrMesh.cpp:964
void parseGridGeneration()
Parse the grid generation algorithm.
Definition CD_AmrMesh.cpp:2654
void interpToNewGrids(MFAMRCellData &a_newData, const MFAMRCellData &a_oldData, const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel, const EBCoarseToFineInterp::Type a_type)
Interpolate data to new grids.
Definition CD_AmrMesh.cpp:2077
void transferCoveredParticlesVoxels(ParticleContainer< P > &a_particlesFrom, ParticleContainer< P > &a_particlesTo, const phase::which_phase &a_phase) const
Function which transferse particles from one particle container to another if they fall inside the EB...
Definition CD_AmrMeshImplem.H:832
void parseMaxEbisBoxSize()
Parse the maximum permitted box size.
Definition CD_AmrMesh.cpp:2616
void alias(Vector< T * > &a_alias, const Vector< RefCountedPtr< T > > &a_data) const
Turn smart-pointer data structure into regular-pointer data structure.
Definition CD_AmrMeshImplem.H:209
EBAMRParticleMesh & getParticleMesh(const std::string a_realm, const phase::which_phase a_phase) const
Get EBAMRParticleMesh operator.
Definition CD_AmrMesh.cpp:3324
bool m_hasGrids
Has grids or not.
Definition CD_AmrMesh.H:2021
Vector< RefCountedPtr< LayoutData< VoFIterator > > > & getVofIterator(const std::string a_realm, const phase::which_phase a_phase) const
Get vof iterators for a Realm and phase. This has the capability of iterating through cut-cells.
Definition CD_AmrMesh.cpp:3192
void deallocate(Vector< T * > &a_data) const
Deallocate data.
Definition CD_AmrMeshImplem.H:161
int getMaxAmrDepth() const
Get maximum permitted amr depth.
Definition CD_AmrMesh.cpp:2964
const Vector< int > & getRefinementRatios() const
Get refinement ratios.
Definition CD_AmrMesh.cpp:3115
void regridOperators(const int a_lmin)
Regrid AMR operators. This is done for all realms.
Definition CD_AmrMesh.cpp:1097
int getRedistributionRadius() const
Get default redistribution radius.
Definition CD_AmrMesh.cpp:3008
const Vector< RefCountedPtr< LevelTiles > > & getLevelTiles(const std::string a_realm) const
Get the tiled space representation.
Definition CD_AmrMesh.cpp:3241
BoxSorting getBoxSorting() const
Get box sorting method.
Definition CD_AmrMesh.cpp:3678
int m_maxEbisBoxSize
Maximum box size for EBIS generation.
Definition CD_AmrMesh.H:1971
int getNumberOfEbGhostCells() const
Get number of ghost cells used for EB grid generation.
Definition CD_AmrMesh.cpp:2997
const EBAMRFAB & getLevelset(const std::string a_realm, const phase::which_phase a_phase) const
Get levelset function, allocated over a grid for a Realm and phase.
Definition CD_AmrMesh.cpp:3308
void regridRealm(const std::string a_realm, const Vector< Vector< int > > &a_procs, const Vector< Vector< Box > > &a_boxes, const int a_lmin)
Regrid a realm. This generates the grids for the realm, but does not do the operators on the realm.
Definition CD_AmrMesh.cpp:3604
void setBaseImplicitFunction(const phase::which_phase a_phase, const RefCountedPtr< BaseIF > &a_baseIF)
Set implicit function for a specific phase. Need e.g. for level-sets.
Definition CD_AmrMesh.cpp:901
int getFinestLevel() const
Get finest grid level.
Definition CD_AmrMesh.cpp:2953
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGridCoFi(const std::string a_realm, const phase::which_phase a_phase) const
Get the coarsened fine EBLevelGrid for a Realm and phase.
Definition CD_AmrMesh.cpp:3274
int m_redistributionRadius
Redistribution radius.
Definition CD_AmrMesh.H:2016
const Vector< ProblemDomain > & getDomains() const
Get domains.
Definition CD_AmrMesh.cpp:3137
int m_numLsfGhostCells
Number of ghost cells to use when writing level-set to grid.
Definition CD_AmrMesh.H:1996
Type
Supported interpolation types.
Definition CD_CellCentroidInterpolation.H:43
Class for handling particle-mesh operations with AMR.
Definition CD_EBAMRParticleMesh.H:52
class for handling surface deposition of particles with EB and AMR.
Definition CD_EBAMRSurfaceDeposition.H:29
Type
Supported interpolation types.
Definition CD_EBCentroidInterpolation.H:43
Type
Type of interpolation methods supported. PWC = Piecewise constant, ignoring the embedded boundary....
Definition CD_EBCoarseToFineInterp.H:42
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25