Mushy Layer
1.0
|
this class manages the various forms of the projection More...
#include <Projector.H>
Public Member Functions | |
Projector () | |
default constructor | |
Projector (const DisjointBoxLayout &a_grids, const DisjointBoxLayout *a_crseGridsPtr, const Box &a_domain, const Real a_dx, Projector *a_finerProj, Projector *a_crseProj, int a_nRefCrse, int a_level, PhysBCUtil &a_physBC) | |
full constructor | |
Projector (const DisjointBoxLayout &a_grids, const DisjointBoxLayout *a_crseGridsPtr, const ProblemDomain &a_domain, const Real a_dx, Projector *a_finerProj, Projector *a_crseProj, int a_nRefCrse, int a_level, PhysBCUtil &a_physBC) | |
full constructor | |
void | copyPressure (Projector &a_proj) |
Copy pressure from another projection object. | |
void | rescalePressure (Real a_oldDt, Real a_newDt) |
May need to rescale pressures after changing the timestep? | |
~Projector () | |
destructor | |
void | define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout *a_crseGridsPtr, const Box &a_domain, const Real a_dx, Projector *a_finerProj, Projector *a_crseProj, int a_nRefCrse, int a_level, PhysBCUtil &a_physBC, bool a_usePiAdvectionBCs=true) |
define function | |
void | define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout *a_crseGridsPtr, const ProblemDomain &a_domain, const Real a_dx, Projector *a_finerProj, Projector *a_crseProj, int a_nRefCrse, int a_level, PhysBCUtil &a_physBC, bool a_usePiAdvectionBCs=true) |
define function | |
void | init (const Projector &a_oldProj) |
initialize new projection with data from old projection | |
void | variableSetUp () |
define static parts | |
void | setCrseProj (Projector *a_crseProj, int nRefCrse) |
set coarse projection operator | |
void | setFineProj (Projector *a_fineProj) |
set fine projection operator | |
void | scalePwithPorosity (bool a_scalePwithPorosity) |
set whether to scale pressure with porosity | |
void | limitSolverCoarsening (bool a_limitSolverCoarsening) |
set solver parameter | |
void | writeCheckpointHeader (HDF5Handle &a_handle) const |
write checkpoint header | |
void | writeCheckpointLevel (HDF5Handle &a_handle) const |
write this class to a checkpoint file for later restart | |
void | readCheckpointHeader (HDF5Handle &a_handle) |
read the checkpoint header | |
void | readCheckpointLevel (HDF5Handle &a_handle) |
read this class from a checkpoint file | |
int | levelMacProject (LevelData< FluxBox > &uEdge, Real a_oldTime, Real a_dt, const RefCountedPtr< LevelData< FArrayBox > > a_porosityPtr, const RefCountedPtr< LevelData< FArrayBox > > a_crsePorosityPtr, const RefCountedPtr< LevelData< FluxBox > > a_porosityEdgePtr, const RefCountedPtr< LevelData< FluxBox > > a_crsePorosityEdgePtr) |
do level MAC projection, correct uEdge. | |
Real | getPhiScale (Real a_dt) |
Scale for . More... | |
void | checkDivergence (LevelData< FluxBox > &a_uEdge) |
Check if the divergence of the velocity is sufficiently small. More... | |
Real | getScale (Real a_scale, Real a_dt) |
Return the scale for . | |
void | setNonSubcycledMACBCs () |
Set up for nonsubcycled solves. | |
void | setSubcycledMACBCs () |
Set up for subcycled solves. | |
int | levelMacProject (LevelData< FluxBox > &a_uEdge, Real a_dt, const RefCountedPtr< LevelData< FArrayBox > > a_crsePressurePtr, const RefCountedPtr< LevelData< FArrayBox > > a_porosityPtr, const RefCountedPtr< LevelData< FArrayBox > > a_crsePorosityPtr, const RefCountedPtr< LevelData< FluxBox > > a_porosityEdgePtr, const RefCountedPtr< LevelData< FluxBox > > a_crsePorosityEdgePtr, bool alreadyhasPhi=false, Real correctScale=1.0) |
Projection function that let's user specify pressure coarse fine boundary condition. | |
void | getCFBC (LevelData< FArrayBox > &crseVelPtr, LevelData< FArrayBox > *a_crseVelPtr, LevelData< FArrayBox > *a_crsePorosityPtr, Real a_dt) |
Get the boundary condition on internal (coarse-fine) boundaries. | |
void | setPressureScalePtr (RefCountedPtr< LevelData< FArrayBox > > a_pressureScalePtr) |
Set pointer to the pressure scaling factor (cell centred) | |
void | setPressureScaleEdgePtr (RefCountedPtr< LevelData< FluxBox > > a_pressureScaleEdgePtr) |
Set pointer to the pressure scaling factor (face centred) | |
void | scaleRHS (LevelData< FArrayBox > &a_rhs, Real a_scale) |
Refactored to use in implicit MAC solve. | |
void | LevelProject (LevelData< FArrayBox > &a_velocity, LevelData< FArrayBox > *a_crseVelPtr, const Real a_newTime, const Real a_dt, const RefCountedPtr< LevelData< FArrayBox > > a_porosityPtr, const RefCountedPtr< LevelData< FArrayBox > > a_crsePorosityPtr, const RefCountedPtr< LevelData< FluxBox > > a_porosityEdgePtr, const RefCountedPtr< LevelData< FluxBox > > a_crsePorosityEdgePtr, const bool a_isViscous) |
do level projection and correct (cell-centered) velocities; | |
void | AdditionalLevelProject (LevelData< FArrayBox > &a_velocity, LevelData< FArrayBox > *a_crseVelPtr, const Real a_newTime, const Real a_dt, const bool a_isViscous) |
For doing a further (constant coefficient) projection. | |
void | doSyncOperations (Vector< LevelData< FArrayBox > * > &a_velocity, Vector< LevelData< FArrayBox > * > &a_lambda, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosityFace, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, const Real a_newTime, const Real a_dtSync) |
do synchronisation operations More... | |
void | initialLevelProject (LevelData< FArrayBox > &a_velocity, LevelData< FArrayBox > *a_crseVelPtr, const Real a_oldTime, const Real a_newTime, const RefCountedPtr< LevelData< FArrayBox > > a_porosityPtr, const RefCountedPtr< LevelData< FArrayBox > > a_crsePorosityPtr, const RefCountedPtr< LevelData< FluxBox > > a_porosityEdgePtr, const RefCountedPtr< LevelData< FluxBox > > a_crsePorosityEdgePtr) |
Initialization functions. More... | |
void | doInitialSyncOperations (Vector< LevelData< FArrayBox > * > &a_vel, Vector< LevelData< FArrayBox > * > &a_lambda, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosityFace, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, const Real a_newTime, const Real a_dtSync) |
void | initialVelocityProject (Vector< LevelData< FArrayBox > * > &a_velocity, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosityFace, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, bool a_homogeneousCFBC=true) |
performs multilevel projection on velocity, modifies velocity in place More... | |
void | initialVelocityProject (Vector< LevelData< FArrayBox > * > &a_velocity, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, AMRMultiGrid< LevelData< FArrayBox > > &a_solver, bool a_homogeneousCFBC=true) |
same as other initialVelocityProject, but uses a pre-defined AMRSolver | |
void | doPostRegridOps (Vector< LevelData< FArrayBox > * > &a_lambda, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosity, const Real a_dt, const Real a_time, const Real a_etaScale=1.0) |
performs the necessary re-initialization after regridding | |
int | getLevel () const |
access functions More... | |
const ProblemDomain & | dProblem () const |
get the problem domain | |
const DisjointBoxLayout & | getBoxes () const |
get the DisjointBoxLayout for this projection operator | |
int | nRefCrse () const |
get the refinement ratio to the coarser level | |
Real | dx () const |
get the cell spacing | |
Real | etaLambda () const |
returns coefficient for volume-discrepancy solve. | |
void | etaLambda (Real a_eta) |
We may want to change the value of during a simulation. | |
Projector * | fineProjPtr () const |
pointer to fine projection operator | |
Projector * | crseProjPtr () const |
pointer to coarse projection operator | |
LevelData< FArrayBox > & | phi () |
returns MAC correction | |
void | getPhi (LevelData< FArrayBox > &a_phi) |
Get . | |
const LevelData< FArrayBox > & | phi () const |
const version of accessor | |
void | gradPhi (LevelData< FArrayBox > &a_gradPhi, int a_dir) const |
returns edge-centered grad(phi) in direction dir | |
void | gradPhi (LevelData< FluxBox > &a_gradPhi) const |
returns all components of grad(phi) (gradPhi should be correct size) | |
LevelData< FArrayBox > & | Pi () |
returns level-projection pressure | |
LevelData< FArrayBox > & | CCrhs () |
RHS for projection, i.e. (cell centred) | |
LevelData< FArrayBox > & | MACrhs () |
RHS for projection solve, i.e. (face centred) | |
LevelData< FArrayBox > & | CCcorrection () |
Correction to make velocity incompressible (cell centred) | |
LevelData< FluxBox > & | MACcorrection () |
Correction to make velocity incompressible (face centred) | |
void | unscaledPi (LevelData< FArrayBox > &unscaledPi, Real a_dt) |
Pressure, not scaled by the timestep. | |
void | gradPi (LevelData< FArrayBox > &a_gradPi, int a_dir) const |
returns grad(pi) in direction dir | |
void | gradPi (LevelData< FArrayBox > &a_gradPi) const |
returns grad(pi) in all directions into SpaceDim-dimensioned gradPi | |
void | gradPiBCs (LevelData< FArrayBox > &a_gradPi, bool extrapBCs=false, bool a_usePhi=false) |
Applys BCs. | |
LevelData< FArrayBox > & | eSync () |
returns synchronization correction | |
const LevelData< FArrayBox > & | eSync () const |
returns synchronization correction (const version) | |
void | grad_eSync (LevelData< FArrayBox > &a_grad_eSync, int a_dir) const |
returns cell-centered grad(eSync) (composite gradient) | |
void | grad_eSync (LevelData< FArrayBox > &a_grad_eSync) const |
returns cell-centered G^{comp}(eSync) | |
LevelData< FArrayBox > & | eLambda () |
returns volume-discrepancy correction | |
const LevelData< FArrayBox > & | eLambda () const |
returns volume-discrepancy correction (const version) | |
void | grad_eLambda (LevelData< FArrayBox > &a_grad_eLambda, int a_dir) const |
returns edge-centered grad(eLambda) in direction dir | |
void | applyFreestreamCorrection (LevelData< FluxBox > &a_advVel, Real scale=1.0) |
correct face centred velocities with the already calculated grad(e_lambda) | |
void | removeFreestreamCorrection (LevelData< FluxBox > &a_advVel) |
Remove freestream correction from face centred velocities. | |
LevelData< FluxBox > & | grad_eLambda () |
returns edge-centered grad(eLambda) in all directions More... | |
const LevelData< FluxBox > & | grad_eLambda () const |
returns edge-centered grad(eLambda) in all directions | |
bool | doSyncProjection () const |
should we do sync projection? | |
bool | doMacSync () const |
should we do volume discrepancy correction? | |
bool | isInitialized () const |
has this object been completely initialized? | |
bool | doQuadInterp () const |
use quadratic interpolation (instead of extrap) for c/f bc | |
QuadCFInterp & | quadCFInterpolator () |
returns predefined quadratic coarse-fine interpolation object | |
const LayoutData< IntVectSet > & | getGridIVS () |
returns predefined intvectset which shows coverage | |
bool | isFinestLevel () const |
is this the finest level? | |
void | isFinestLevel (bool a_finest_level) |
set whether this is the finest level | |
void | verbosity (int a_verbosity) |
set verbosity | |
int | verbosity () const |
Set fine proj ptr. More... | |
void | setPhysBC (PhysBCUtil &a_bc) |
sets physBCs | |
PhysBCUtil * | getPhysBCPtr () const |
pointer to the boundary condition object | |
void | computeVDCorrection (Vector< LevelData< FArrayBox > * > &a_lambda, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosity, const Real a_newTime, const Real a_dtSync) |
Version which defines it's own amr multigrid. | |
void | applyMacCorrection (LevelData< FluxBox > &a_uEdge, LevelData< FArrayBox > *crseBCDataPtr, Real phiScale=1.0) |
returns uEdge - G(phi_mac); C/F BC is CFscale*Pi(crse) More... | |
Public Attributes | |
Real | m_sumVDrhs |
Sum of RHS for VD correction. | |
Protected Member Functions | |
void | applyCCcorrection (LevelData< FArrayBox > &a_velocity, LevelData< FArrayBox > *a_pressure, LevelData< FArrayBox > *a_pressureScalePtr, const Real scale) |
vel := vel - scale*G_CC(pi) More... | |
void | doSyncProjection (Vector< LevelData< FArrayBox > * > &a_velocity, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, const Real a_newTime, const Real a_dtSync, AMRMultiGrid< LevelData< FArrayBox > > &a_solver) |
do sync projection with this level as lBase; corrects velocity in place | |
void | initialSyncProjection (Vector< LevelData< FArrayBox > * > &a_velocity, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, const Real a_newTime, const Real a_dtSync, AMRMultiGrid< LevelData< FArrayBox > > &a_solver) |
perform initial sync projection | |
void | computeVDCorrection (Vector< LevelData< FArrayBox > * > &a_lambda, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosity, const Real a_newTime, const Real a_dtSync, AMRMultiGrid< LevelData< FArrayBox > > &a_solver) |
do volume discrepancy solve, with this level as lBase | |
void | applySyncCorrection (Vector< LevelData< FArrayBox > * > &a_velocity, Vector< RefCountedPtr< LevelData< FArrayBox > > > &a_porosity, const Real scale, LevelData< FArrayBox > *crseCorr) |
Compute volume discrepancy with this level as lMax. More... | |
void | computeGrad_eLambda (Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosity) |
compute composite gradient of eLambda – does this recursively | |
void | defineMultiGrid (AMRMultiGrid< LevelData< FArrayBox > > &a_solver, const Vector< LevelData< FArrayBox > * > &a_vel, Vector< RefCountedPtr< LevelData< FluxBox > > > &a_porosity, bool a_freestreamSolver) |
rescales composite gradient of eLambda, recursively all finer levels More... | |
void | defineSolverMGlevel (const DisjointBoxLayout &a_grids, const DisjointBoxLayout *a_crseGridsPtr) |
Define multigrid solver for freestream preservation solve. More... | |
void | setSolverParameters () |
Set things like tolerance, number of smoothing steps etc. | |
void | makeBottomSolvers () |
Setup bototm solvers. | |
void | defineSolverMGLevelVCOp (const RefCountedPtr< LevelData< FArrayBox > > a_porosityPtr, const RefCountedPtr< LevelData< FArrayBox > > a_crsePorosityPtr, const RefCountedPtr< LevelData< FluxBox > > a_porosityEdgePtr, const RefCountedPtr< LevelData< FluxBox > > a_crsePorosityEdgePtr, bool cellCentred=false, Real beta=-1) |
define multigrid solver for this level More... | |
int | solveMGlevel (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > *a_phiCoarsePtr, const LevelData< FArrayBox > &a_rhs, const RefCountedPtr< LevelData< FArrayBox > > a_porosityPtr, const RefCountedPtr< LevelData< FArrayBox > > a_crsePressureScalePtr, const RefCountedPtr< LevelData< FluxBox > > a_porosityEdgePtr, const RefCountedPtr< LevelData< FluxBox > > a_crsePorosityEdgePtr, bool cellCentred=false) |
solve for pressure on this level | |
void | solveMGlevel (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > *a_phiCoarsePtr, const LevelData< FArrayBox > &a_rhs, bool cellCentred) |
For when we don't have any pressure scale. | |
void | postRestart () |
takes care of setting ghost cell values after restarting | |
this class manages the various forms of the projection
this class performs the various cell-centered projections required by the IAMR algorithm. These include the single-level operators levelMac and levelProject, as well as the multilevel sync projection and volume discrepancy solves.
|
protected |
vel := vel - scale*G_CC(pi)
where is the scale
void Projector::applyMacCorrection | ( | LevelData< FluxBox > & | a_uEdge, |
LevelData< FArrayBox > * | crseBCDataPtr, | ||
Real | phiScale = 1.0 |
||
) |
returns uEdge - G(phi_mac); C/F BC is CFscale*Pi(crse)
|
protected |
Compute volume discrepancy with this level as lMax.
apply sync correction vel = vel - scale*G^{comp} eSync (recursively applies correction to all finer levels)
Check if the divergence of the velocity is sufficiently small.
If not, increase the number of smoothing steps
|
protected |
rescales composite gradient of eLambda, recursively all finer levels
rescales grad(eLambda) by (a_dx_lbase)/m_dx for this level and recursively calls this function for finer levels. This is done for stability reasons.defines Multilevel AMRSolver used for composite projections uses grids in vel, lBase is this level. a_freestream solver is true if solver will be used for freestream preservation solve.
|
protected |
Define multigrid solver for freestream preservation solve.
defines multigrid solver for this level
|
protected |
define multigrid solver for this level
variable coefficient version to deal with variable porosity
void Projector::doInitialSyncOperations | ( | Vector< LevelData< FArrayBox > * > & | a_vel, |
Vector< LevelData< FArrayBox > * > & | a_lambda, | ||
Vector< RefCountedPtr< LevelData< FluxBox > > > & | a_porosityFace, | ||
Vector< RefCountedPtr< LevelData< FArrayBox > > > & | a_porosity, | ||
const Real | a_newTime, | ||
const Real | a_dtSync | ||
) |
defines AMRSolver and calls initial sync projection and freestream preservation solves
void Projector::doSyncOperations | ( | Vector< LevelData< FArrayBox > * > & | a_velocity, |
Vector< LevelData< FArrayBox > * > & | a_lambda, | ||
Vector< RefCountedPtr< LevelData< FluxBox > > > & | a_porosityFace, | ||
Vector< RefCountedPtr< LevelData< FArrayBox > > > & | a_porosity, | ||
const Real | a_newTime, | ||
const Real | a_dtSync | ||
) |
do synchronisation operations
defines AMRSolver and calls sync projection and freestream preservation solve. assumes all physical and copy BC's already set (does C/F BC's however)
int Projector::getLevel | ( | ) | const |
access functions
first define quick-n-easy access functions
get the level this Projector is defined for
Scale for .
Important to get this right for coarse-fine interface matching
returns edge-centered grad(eLambda) in all directions
non-constant reference returned because of need to rescale.
void Projector::initialLevelProject | ( | LevelData< FArrayBox > & | a_velocity, |
LevelData< FArrayBox > * | a_crseVelPtr, | ||
const Real | a_oldTime, | ||
const Real | a_newTime, | ||
const RefCountedPtr< LevelData< FArrayBox > > | a_porosityPtr, | ||
const RefCountedPtr< LevelData< FArrayBox > > | a_crsePorosityPtr, | ||
const RefCountedPtr< LevelData< FluxBox > > | a_porosityEdgePtr, | ||
const RefCountedPtr< LevelData< FluxBox > > | a_crsePorosityEdgePtr | ||
) |
Initialization functions.
performs initial level projection – phys and copy BC's should be preset
void Projector::initialVelocityProject | ( | Vector< LevelData< FArrayBox > * > & | a_velocity, |
Vector< RefCountedPtr< LevelData< FluxBox > > > & | a_porosityFace, | ||
Vector< RefCountedPtr< LevelData< FArrayBox > > > & | a_porosity, | ||
bool | a_homogeneousCFBC = true |
||
) |
performs multilevel projection on velocity, modifies velocity in place
dtSync is used as a scaling factor for coarse-level BC's, not used in actual projection. if homogeneousCFBC is false, C/F bc for l_base is taken from coarse-level eSync, o/w use l_base-1 correction field = 0. if applyCorrection == false, then correction is computed, but not applied (can be useful as a diagnostic). physical and copy BC's should already be set before this is called.
int Projector::verbosity | ( | ) | const |
Set fine proj ptr.
returns verbosity