BISICLES AMR ice sheet model
0.9
|
Basic Sigma fourth-order coordinate system on an AMR level. More...
#include <LevelSigmaCS.H>
Public Types | |
enum | ThicknessInterpolationMethod { STD_THICKNESS_INTERPOLATION_METHOD = 0, SMOOTH_SURFACE_THICKNESS_INTERPOLATION_METHOD = 1, MAX_THICKNESS_INTERPOLATION_METHOD = 2 } |
choice of thickness interpolation method in interpFromCoarse More... | |
Public Member Functions | |
LevelSigmaCS () | |
default constructor More... | |
LevelSigmaCS (const DisjointBoxLayout &a_grids, const RealVect &a_dx, const IntVect &a_ghostVect=IntVect::Unit) | |
defining constructor – calls matching define function More... | |
LevelSigmaCS (const DisjointBoxLayout &a_grids, const RealVect &a_dx, const LevelSigmaCS &a_fineCS, int a_nRef) | |
define as a coarsening of fineCS More... | |
virtual | ~LevelSigmaCS () |
void | define (const DisjointBoxLayout &a_grids, const RealVect &a_dx, const IntVect &a_ghostVect=IntVect::Unit) |
void | define (const LevelSigmaCS &a_fineCS, int a_nRef) |
define as a coarsening of fineCS by nRef More... | |
RealVect | realCoord (const RealVect &a_Xi, const DataIndex &a_index) const |
RealVect | mappedCoord (const RealVect &a_x, const DataIndex &a_index) const |
given coordinate in real space, return location in mapped space More... | |
void | getNodeRealCoordinates (LevelData< NodeFArrayBox > &a_nodeCoords) const |
return Cartesian XYZ locations of nodes More... | |
RealVect | dx () const |
LevelData< FArrayBox > & | getH () |
returns modifiable cell-centered H (ice sheet thickness) for this More... | |
const LevelData< FArrayBox > & | getH () const |
returns const reference to cell-centered H (ice sheet thickness) More... | |
LevelData< FluxBox > & | getFaceH () |
returns modifiable face-centered H More... | |
const LevelData< FluxBox > & | getFaceH () const |
returns const reference to face-centered H More... | |
const LevelData< FArrayBox > & | getTopography () const |
returns a const-reference to the cell-centered topography More... | |
LevelData< FArrayBox > & | getTopography () |
returns a modifiable reference to the cell-centered topography More... | |
void | setTopography (const LevelData< FArrayBox > &a_topography) |
sets the base height. More... | |
void | setBackgroundSlope (const RealVect &a_backgroundSlope) |
set the unit topography difference More... | |
const RealVect & | getBackgroundSlope () const |
void | setSurfaceHeight (const LevelData< FArrayBox > &a_surface) |
sets the cell-centered surface height More... | |
void | setGradSurface (const LevelData< FArrayBox > &a_gradSurface) |
sets the cell-centered surface gradient More... | |
void | setGradSurfaceFace (const LevelData< FluxBox > &a_gradSurface) |
sets the face-centred surface gradient More... | |
const LevelData< FArrayBox > & | getSurfaceHeight () const |
returns the cell-centred surface elevation More... | |
const LevelData< FArrayBox > & | getGradSurface () const |
returns the cell-centred surface gradient More... | |
const LevelData< FluxBox > & | getGradSurfaceFace () const |
returns the face-centred surface gradient More... | |
const LevelData< FArrayBox > & | getThicknessOverFlotation () const |
void | getSurfaceHeight (LevelData< FArrayBox > &a_zSurface) const |
computes ice surface height More... | |
const LevelData< BaseFab< int > > & | getFloatingMask () const |
returns floating mask specifying whether ice is grounded or floating More... | |
void | recomputeGeometry (const LevelSigmaCS *a_crseCoords, const int a_refRatio) |
recomputes quantities which are dependent on cell-centered H. More... | |
void | recomputeGeometryFace (const LevelSigmaCS *a_crseCoords, const int a_refRatio) |
recomputes quantities which are dependent on face-centered H. More... | |
const LevelData< FArrayBox > & | deltaFactors () const |
return cell-centered ![]() ![]() | |
const LevelData< FluxBox > & | faceDeltaFactors () const |
return face-centered ![]() ![]() | |
Real | seaLevel () const |
returns sea level More... | |
void | setSeaLevel (const Real &a_seaLevel) |
set sea level More... | |
Real | iceDensity () const |
ice density – this is probably an interim approach More... | |
Real | waterDensity () const |
seawater density – probably an interim approach More... | |
Real | gravity () const |
acceleration due to gravity – probably an interim approach More... | |
void | setIceDensity (const Real &a_iceDensity) |
interim or not, we also need to set these More... | |
void | setWaterDensity (const Real &a_waterDensity) |
void | setGravity (const Real &a_gravity) |
bool | isDefined () const |
returns true if this object has been defined. More... | |
LevelSigmaCS * | makeCoarser (int a_coarsenFactor) const |
const DisjointBoxLayout & | grids () const |
access function More... | |
const LayoutData< bool > & | anyFloating () const |
void | interpFromCoarse (const LevelSigmaCS &a_crseCoords, const int a_refinementRatio, const bool a_interpolateTopography=true, const bool a_interpolateThickness=true, const bool a_preserveMask=true, const bool a_interpolateTopographyGhost=true, const bool a_interpolateThicknessGhost=true, const bool a_preserveMaskGhost=true, int a_thicknessInterpolationMethod=STD_THICKNESS_INTERPOLATION_METHOD) |
fill the data in this by interpolation from a coarser level More... | |
void | exchangeTopography () |
void | unitShiftExchange (LevelData< FArrayBox > &a_level) |
IntVect | ghostVect () const |
accessor for ghost vector More... | |
const Vector< Real > & | getFaceSigma () const |
const Vector< Real > & | getDSigma () const |
const Vector< Real > & | getSigma () const |
void | setFaceSigma (const Vector< Real > &a_faceSigma) |
Protected Member Functions | |
void | setDefaultValues () |
void | computeDeltaFactors () |
(re)compute delta factors More... | |
void | computeSurface (const LevelSigmaCS *a_crseCoords, const int a_refRatio) |
(re)compute surface elevation and gradients More... | |
void | computeFloatingMask (const LevelData< FArrayBox > &a_surface) |
sets floating mask based on surface height More... | |
Protected Attributes | |
DisjointBoxLayout | m_grids |
DisjointBoxLayout | m_Hgrids |
"flattened" grids (in 3d) for quantities like H More... | |
RealVect | m_dx |
IntVect | m_ghostVect |
ghost vector More... | |
RealVect | m_unitShift |
RealVect | m_backgroundSlope |
LevelData< FArrayBox > | m_topography |
cell-centered topography More... | |
LevelData< FArrayBox > | m_H |
cell-centered ice thickness More... | |
LevelData< FluxBox > | m_faceH |
face-centered ice thickness More... | |
LevelData< FArrayBox > | m_surface |
cell-centered surface elevation More... | |
LevelData< FArrayBox > | m_gradSurface |
cell-centered gradient of surface elevation More... | |
LevelData< FluxBox > | m_gradSurfaceFace |
face-centered gradient of surface elevation More... | |
LevelData< FArrayBox > | m_thicknessOverFlotation |
cell-centered thickness over flotation More... | |
LevelData< BaseFab< int > > | m_floatingMask |
floating mask More... | |
LayoutData< bool > | m_anyFloating |
is anything floating? More... | |
LevelData< FArrayBox > | m_deltaFactors |
LevelData< FluxBox > | m_faceDeltaFactors |
bool | m_isDefined |
has this object been defined? More... | |
Real | m_seaLevel |
Real | m_iceDensity |
ice density More... | |
Real | m_waterDensity |
seawater density – probably an interim approach More... | |
Real | m_gravity |
acceleration due to gravity – probably an interim approach More... | |
Vector< Real > | m_faceSigma |
Vector< Real > | m_sigma |
Vector< Real > | m_dSigma |
Basic Sigma fourth-order coordinate system on an AMR level.
The SigmaCS class implements Sigma-coordinates for ice sheets on a DisjointBoxLayout
choice of thickness interpolation method in interpFromCoarse
The original (default) method, STD_THICKNESS_INTERPOLATION_METHOD, interpolates the thickness directly, which can lead to a smooth thickness on a rough bedrock, and henec a rough surface
An alternative method, SMOOTH_SURFACE_THICKNESS_INTERPOLATION_METHOD, interpolates the upper and lower surfaces, but replaces the lower surface on grounded (on the coarse mesh) cells with the bedrock data. This will be conservative iff the relationship between fine and coarse bedrock data is conservative
Enumerator | |
---|---|
STD_THICKNESS_INTERPOLATION_METHOD | |
SMOOTH_SURFACE_THICKNESS_INTERPOLATION_METHOD | |
MAX_THICKNESS_INTERPOLATION_METHOD |
LevelSigmaCS::LevelSigmaCS | ( | ) |
LevelSigmaCS::LevelSigmaCS | ( | const DisjointBoxLayout & | a_grids, |
const RealVect & | a_dx, | ||
const IntVect & | a_ghostVect = IntVect::Unit |
||
) |
defining constructor – calls matching define function
defining constructor
References define(), and setDefaultValues().
LevelSigmaCS::LevelSigmaCS | ( | const DisjointBoxLayout & | a_grids, |
const RealVect & | a_dx, | ||
const LevelSigmaCS & | a_fineCS, | ||
int | a_nRef | ||
) |
define as a coarsening of fineCS
cell-centered topography
cell-centered gradient of surface elevation
cell-centered ice thickness
cell-centered surface elevation
cell-centered
cell-centered
face-centered ice thickness
face-centered gradient of surface elevation
face-centered
References averageAllDim(), define(), ghostVect(), horizontalAverage(), m_deltaFactors, m_dSigma, m_faceDeltaFactors, m_faceH, m_faceSigma, m_gradSurface, m_gradSurfaceFace, m_grids, m_H, m_seaLevel, m_sigma, m_surface, m_thicknessOverFlotation, m_topography, seaLevel(), and setDefaultValues().
|
virtual |
Destructor.
|
inline |
References m_anyFloating.
Referenced by IceUtility::defineRHS(), and IceUtility::setFloatingBasalFriction().
|
protected |
(re)compute delta factors
References m_deltaFactors, m_dx, m_grids, m_H, and m_topography.
Referenced by recomputeGeometry(), and recomputeGeometryFace().
|
protected |
sets floating mask based on surface height
References m_anyFloating, m_backgroundSlope, m_floatingMask, m_grids, m_H, m_iceDensity, m_seaLevel, m_topography, and m_waterDensity.
Referenced by computeSurface().
|
protected |
(re)compute surface elevation and gradients
References computeFloatingMask(), dx(), getSurfaceHeight(), iceDensity(), m_backgroundSlope, m_dx, m_floatingMask, m_gradSurface, m_gradSurfaceFace, m_grids, m_H, m_seaLevel, m_surface, m_thicknessOverFlotation, m_topography, column_thermodynamics::rhoi, column_thermodynamics::rhoo, seaLevel(), and waterDensity().
Referenced by recomputeGeometry(), and recomputeGeometryFace().
void LevelSigmaCS::define | ( | const DisjointBoxLayout & | a_grids, |
const RealVect & | a_dx, | ||
const IntVect & | a_ghostVect = IntVect::Unit |
||
) |
References m_anyFloating, m_deltaFactors, m_dx, m_faceDeltaFactors, m_faceH, m_floatingMask, m_ghostVect, m_gradSurface, m_gradSurfaceFace, m_grids, m_H, m_Hgrids, m_isDefined, m_surface, m_thicknessOverFlotation, and m_topography.
Referenced by define(), LevelSigmaCS(), and makeCoarser().
void LevelSigmaCS::define | ( | const LevelSigmaCS & | a_fineCS, |
int | a_nRef | ||
) |
define as a coarsening of fineCS by nRef
References averageAllDim(), averageAllDimFace(), define(), ghostVect(), horizontalAverage(), horizontalAverageFace(), m_deltaFactors, m_dSigma, m_dx, m_faceDeltaFactors, m_faceH, m_faceSigma, m_gradSurface, m_gradSurfaceFace, m_grids, m_H, m_seaLevel, m_sigma, m_topography, and seaLevel().
const LevelData< FArrayBox > & LevelSigmaCS::deltaFactors | ( | ) | const |
return cell-centered and
,
return cell-centered ${{x}}$ and ${{y}}$,
stored as components 0 and 1, respectively.
References m_deltaFactors.
Referenced by computeCCDerivatives(), and getFloatingMask().
|
inline |
References getFaceH(), getH(), getTopography(), m_dx, and setTopography().
Referenced by computeCCDerivatives(), computeFCDerivatives(), ConstitutiveRelation::computeStrainRateInvariantFace(), computeSurface(), defineLevelSigmaCS(), HumpIBC::initializeIceGeometry(), VieliPayneIBC::initializeIceGeometry(), MarineIBC::initializeIceGeometry(), LevelDataTemperatureIBC::initializeIceInternalEnergy(), PythonInterface::PythonIceTemperatureIBC::initializeIceInternalEnergy(), interpFromCoarse(), PythonInterface::PythonIBC::modifyFaceVelocity(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), PythonInterface::PythonIBC::modifyVelocityRHS(), VieliPayneIBC::modifyVelocityRHS(), VieliPayneBCFunction::operator()(), MarineIBC::regridIceGeometry(), FortranInterfaceBasalFriction::setBasalFriction(), GaussianBumpFriction::setBasalFriction(), twistyStreamFriction::setBasalFriction(), singularStreamFriction::setBasalFriction(), LevelDataBasalFriction::setBasalFriction(), MultiLevelDataBasalFriction::setBasalFriction(), sinusoidalFriction::setBasalFriction(), PythonInterface::PythonBasalFriction::setBasalFriction(), LevelDataMuCoefficient::setMuCoefficient(), MultiLevelDataMuCoefficient::setMuCoefficient(), PythonInterface::PythonMuCoefficient::setMuCoefficient(), PythonInterface::PythonVelocitySolver::solve(), AmrIce::updateGeometry(), AmrIce::updateInternalEnergy(), VieliPayneIBC::velocityGhostBC(), IceThicknessIBC::velocityGhostBC(), MarineIBC::velocitySolveBC(), and WriteSigmaMappedAMRHierarchyHDF5().
void LevelSigmaCS::exchangeTopography | ( | ) |
References m_topography.
Referenced by defineLevelSigmaCS(), FortranInterfaceIBC::initializeIceGeometry(), and interpFromCoarse().
const LevelData< FluxBox > & LevelSigmaCS::faceDeltaFactors | ( | ) | const |
return face-centered and
return face-centered ${{x}}$ and ${{y}}$
stored as components 0 and 1, respectively.
References m_faceDeltaFactors.
Referenced by computeFCDerivatives(), and getFloatingMask().
|
inline |
References getGradSurface(), getGradSurfaceFace(), getSurfaceHeight(), getThicknessOverFlotation(), m_backgroundSlope, setGradSurface(), setGradSurfaceFace(), and setSurfaceHeight().
Referenced by AmrIce::writeAMRPlotFile().
|
inline |
References m_dSigma.
Referenced by L1L2ConstitutiveRelation::computeFaceMu(), L1L2ConstitutiveRelation::computeMu(), and AmrIce::updateInternalEnergy().
LevelData< FluxBox > & LevelSigmaCS::getFaceH | ( | ) |
returns modifiable face-centered H
returns modifiable cell-centered H (ice sheet thickness) for this
References m_faceH.
Referenced by DamageConstitutiveRelation::computeFaceMu(), L1L2ConstitutiveRelation::computeFaceMu(), FASIceViscouseTensorOp::computeMu(), PetscIceSolver::computeMu(), dx(), VieliPayneBCFunction::operator()(), IceNonlinearViscousTensor::setState(), and AmrIce::tagCellsLevel().
const LevelData< FluxBox > & LevelSigmaCS::getFaceH | ( | ) | const |
returns const reference to face-centered H
returns const reference to face-centered H (ice sheet thickness)
References m_faceH.
|
inline |
|
inline |
returns floating mask specifying whether ice is grounded or floating
following mask values are defined in IceConstants.H floatingmaskval – ice shelf (nonzero thickness floating on water) groundedmaskval – bottom of ice sheet is grounded. openseamaskval – zero ice thickness, base height below sea level openlandmaskvel – zero ice thickness, base height above sea level
References deltaFactors(), faceDeltaFactors(), m_floatingMask, recomputeGeometry(), and recomputeGeometryFace().
Referenced by AmrIce::advectIceFrac(), MaskedCalvingModel::applyCriterion(), CrevasseCalvingModel::applyCriterion(), DeglaciationCalvingModelA::applyCriterion(), DomainEdgeCalvingModel::applyCriterion(), ProximityCalvingModel::applyCriterion(), DeglaciationCalvingModelB::applyCriterion(), DamageCalvingModel::applyCriterion(), ThicknessCalvingModel::applyCriterion(), MaximumExtentCalvingModel::applyCriterion(), FlotationCalvingModel::applyCriterion(), CliffCollapseCalvingModel::applyCriterion(), VariableRateCalvingModel::applyCriterion(), RateProportionalToSpeedCalvingModel::applyCriterion(), BasalFrictionPowerLaw::computeAlpha(), IceUtility::computeC0(), FASIceViscouseTensorOp::computeMu(), IceUtility::defineRHS(), IceUtility::eliminateFastIce(), IceUtility::eliminateRemoteIce(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), InverseVerticallyIntegratedVelocitySolver::nDoF(), AmrIce::restart(), IceUtility::setFloatingBasalFriction(), AmrIce::solveVelocityField(), AmrIce::tagCellsLevel(), AmrIce::updateGeometry(), AmrIce::updateGroundingLineProximity(), AmrIce::updateInternalEnergy(), and AmrIce::writeAMRPlotFile().
const LevelData< FArrayBox > & LevelSigmaCS::getGradSurface | ( | ) | const |
returns the cell-centred surface gradient
References m_gradSurface.
Referenced by L1L2ConstitutiveRelation::computeFaceMu(), IceUtility::defineRHS(), and getBackgroundSlope().
const LevelData< FluxBox > & LevelSigmaCS::getGradSurfaceFace | ( | ) | const |
returns the face-centred surface gradient
References m_gradSurfaceFace.
Referenced by L1L2ConstitutiveRelation::computeFaceMu(), and getBackgroundSlope().
LevelData< FArrayBox > & LevelSigmaCS::getH | ( | ) |
returns modifiable cell-centered H (ice sheet thickness) for this
References m_H.
Referenced by AmrIce::advectIceFrac(), IceUtility::computeA(), PressureLimitedBasalFrictionRelation::computeAlpha(), IceUtility::computeC0(), computeCCDerivatives(), L1L2ConstitutiveRelation::computeFaceMu(), CrevasseCalvingModel::computeStressMeasure(), DomainDiagnosticData::computeTotalIce(), AmrIce::computeTotalIce(), defineLevelSigmaCS(), IceUtility::defineRHS(), dx(), IceUtility::eliminateFastIce(), IceUtility::eliminateRemoteIce(), IceUtility::extrapVelocityToMargin(), AmrIce::getIceThickness(), AmrIce::implicitThicknessCorrection(), AmrIce::incrementIceThickness(), MultiLevelDataIBC::initializeIceGeometry(), HumpIBC::initializeIceGeometry(), LevelDataIBC::initializeIceGeometry(), VieliPayneIBC::initializeIceGeometry(), PythonInterface::PythonIBC::initializeIceGeometry(), FortranInterfaceIBC::initializeIceGeometry(), MarineIBC::initializeIceGeometry(), PythonInterface::PythonIceTemperatureIBC::initializeIceInternalEnergy(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), PythonInterface::PythonIBC::modifyVelocityRHS(), VieliPayneIBC::modifyVelocityRHS(), IceUtility::multiplyByGroundedFraction(), InverseVerticallyIntegratedVelocitySolver::nDoF(), LevelDataIBC::new_thicknessIBC(), AmrIce::readCheckpointFile(), PythonInterface::PythonIBC::regridIceGeometry(), AmrIce::restart(), PythonInterface::PythonBasalFriction::setBasalFriction(), MultiLevelDataIBC::setGeometryBCs(), LevelDataIBC::setGeometryBCs(), PythonInterface::PythonIBC::setGeometryBCs(), FortranInterfaceIBC::setGeometryBCs(), MarineIBC::setGeometryBCs(), PythonInterface::PythonMuCoefficient::setMuCoefficient(), PythonInterface::PythonVelocitySolver::solve(), AmrIce::tagCellsLevel(), AMRDamage::timestep(), AmrIce::timeStep(), AmrIce::updateCoordSysWithNewThickness(), AmrIce::updateGeometry(), AmrIce::updateInternalEnergy(), MarineIBC::velocitySolveBC(), AmrIce::writeAMRPlotFile(), AmrIce::writeCheckpointFile(), and WriteSigmaMappedAMRHierarchyHDF5().
const LevelData< FArrayBox > & LevelSigmaCS::getH | ( | ) | const |
returns const reference to cell-centered H (ice sheet thickness)
References m_H.
void LevelSigmaCS::getNodeRealCoordinates | ( | LevelData< NodeFArrayBox > & | a_nodeCoords | ) | const |
return Cartesian XYZ locations of nodes
nodeCoords should have 3 coordinates
nodeCoords should have dimension returned by dimension()
References m_dx, m_grids, and realCoord().
Referenced by WriteSigmaMappedAMRHierarchyHDF5().
|
inline |
const LevelData< FArrayBox > & LevelSigmaCS::getSurfaceHeight | ( | ) | const |
returns the cell-centred surface elevation
References m_surface.
Referenced by CliffCollapseCalvingModel::applyCriterion(), IceUtility::computeC0(), computeSurface(), IceUtility::defineRHS(), IceUtility::extrapVelocityToMargin(), getBackgroundSlope(), interpFromCoarse(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), PythonInterface::PythonIBC::modifyVelocityRHS(), InverseVerticallyIntegratedVelocitySolver::nDoF(), and AmrIce::writeAMRPlotFile().
void LevelSigmaCS::getSurfaceHeight | ( | LevelData< FArrayBox > & | a_zSurface | ) | const |
computes ice surface height
returns ice surface height
if ice is grounded, then this equals z_base + H if ice is floating, then this equals seaLevel + H*(1 - rho_ice/rho_water)
if ice is grounded, then this equals z_base + H if ice is floating, then this equals H*(1 - rho_ice/rho_water)
References m_grids, m_H, m_iceDensity, m_seaLevel, m_topography, and m_waterDensity.
const LevelData< FArrayBox > & LevelSigmaCS::getThicknessOverFlotation | ( | ) | const |
const LevelData< FArrayBox > & LevelSigmaCS::getTopography | ( | ) | const |
returns a const-reference to the cell-centered topography
References m_topography.
Referenced by AmrIce::advectIceFrac(), CrevasseCalvingModel::applyCriterion(), IceUtility::computeC0(), DamageConstitutiveRelation::computeFaceMu(), defineLevelSigmaCS(), dx(), IceUtility::extrapVelocityToMargin(), MultiLevelDataIBC::initializeIceGeometry(), HumpIBC::initializeIceGeometry(), LevelDataIBC::initializeIceGeometry(), VieliPayneIBC::initializeIceGeometry(), PythonInterface::PythonIBC::initializeIceGeometry(), FortranInterfaceIBC::initializeIceGeometry(), MarineIBC::initializeIceGeometry(), PythonInterface::PythonIceTemperatureIBC::initializeIceInternalEnergy(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), IceUtility::multiplyByGroundedFraction(), InverseVerticallyIntegratedVelocitySolver::nDoF(), LevelDataIBC::new_thicknessIBC(), MultiLevelDataIBC::regridIceGeometry(), LevelDataIBC::regridIceGeometry(), PythonInterface::PythonIBC::regridIceGeometry(), FortranInterfaceIBC::regridIceGeometry(), MarineIBC::regridIceGeometry(), PythonInterface::PythonBasalFriction::setBasalFriction(), MultiLevelDataIBC::setGeometryBCs(), LevelDataIBC::setGeometryBCs(), PythonInterface::PythonIBC::setGeometryBCs(), FortranInterfaceIBC::setGeometryBCs(), MarineIBC::setGeometryBCs(), PythonInterface::PythonMuCoefficient::setMuCoefficient(), PythonInterface::PythonVelocitySolver::solve(), AmrIce::tagCellsLevel(), AmrIce::timeStep(), AmrIce::updateGeometry(), AmrIce::writeAMRPlotFile(), AmrIce::writeCheckpointFile(), and WriteSigmaMappedAMRHierarchyHDF5().
LevelData< FArrayBox > & LevelSigmaCS::getTopography | ( | ) |
returns a modifiable reference to the cell-centered topography
References m_topography.
|
inline |
accessor for ghost vector
References m_ghostVect.
Referenced by define(), MultiLevelDataIBC::initializeIceGeometry(), LevelDataIBC::initializeIceGeometry(), LevelSigmaCS(), and WriteSigmaMappedAMRHierarchyHDF5().
|
inline |
acceleration due to gravity – probably an interim approach
References m_gravity.
Referenced by IceUtility::computeA(), PressureLimitedBasalFrictionRelation::computeAlpha(), DamageConstitutiveRelation::computeFaceMu(), L1L2ConstitutiveRelation::computeFaceMu(), BennCalvingModel::computeRemnant(), VdVCalvingModel::computeRemnant(), IceUtility::defineRHS(), VerticalConductionInternalEnergyIBC::initializeIceInternalEnergy(), PythonInterface::PythonIBC::modifyVelocityRHS(), VieliPayneIBC::modifyVelocityRHS(), VieliPayneBCFunction::operator()(), AmrIce::timeStep(), AmrIce::updateInternalEnergy(), and MarineIBC::velocitySolveBC().
|
inline |
access function
References m_grids.
Referenced by MaskedCalvingModel::applyCriterion(), CrevasseCalvingModel::applyCriterion(), DeglaciationCalvingModelA::applyCriterion(), DomainEdgeCalvingModel::applyCriterion(), ProximityCalvingModel::applyCriterion(), DeglaciationCalvingModelB::applyCriterion(), DamageCalvingModel::applyCriterion(), ThicknessCalvingModel::applyCriterion(), MaximumExtentCalvingModel::applyCriterion(), FlotationCalvingModel::applyCriterion(), CliffCollapseCalvingModel::applyCriterion(), VariableRateCalvingModel::applyCriterion(), RateProportionalToSpeedCalvingModel::applyCriterion(), IceUtility::computeA(), BennCalvingModel::computeRemnant(), VdVCalvingModel::computeRemnant(), CrevasseCalvingModel::computeStressMeasure(), defineLevelSigmaCS(), MultiLevelDataIBC::initializeIceGeometry(), HumpIBC::initializeIceGeometry(), LevelDataIBC::initializeIceGeometry(), PythonInterface::PythonIBC::initializeIceGeometry(), FortranInterfaceIBC::initializeIceGeometry(), MarineIBC::initializeIceGeometry(), VerticalConductionInternalEnergyIBC::initializeIceInternalEnergy(), LevelDataTemperatureIBC::initializeIceInternalEnergy(), PythonInterface::PythonIBC::modifyFaceVelocity(), PythonInterface::PythonIBC::modifyVelocityRHS(), LevelDataIBC::new_thicknessIBC(), PythonInterface::PythonIBC::regridIceGeometry(), MarineIBC::regridIceGeometry(), ConstantIceTemperatureIBC::setIceInternalEnergyBC(), ReflectionIceInternalEnergyIBC::setIceInternalEnergyBC(), PythonInterface::PythonIceTemperatureIBC::setIceInternalEnergyBC(), and MarineIBC::velocitySolveBC().
|
inline |
ice density – this is probably an interim approach
References m_iceDensity.
Referenced by CrevasseCalvingModel::applyCriterion(), IceUtility::computeA(), PressureLimitedBasalFrictionRelation::computeAlpha(), DamageConstitutiveRelation::computeFaceMu(), L1L2ConstitutiveRelation::computeFaceMu(), BennCalvingModel::computeRemnant(), VdVCalvingModel::computeRemnant(), computeSurface(), IceUtility::defineRHS(), VerticalConductionInternalEnergyIBC::initializeIceInternalEnergy(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), PythonInterface::PythonIBC::modifyVelocityRHS(), VieliPayneIBC::modifyVelocityRHS(), IceUtility::multiplyByGroundedFraction(), LevelDataIBC::new_thicknessIBC(), VieliPayneBCFunction::operator()(), AmrIce::restart(), AmrIce::timeStep(), AmrIce::updateInternalEnergy(), and MarineIBC::velocitySolveBC().
void LevelSigmaCS::interpFromCoarse | ( | const LevelSigmaCS & | a_crseCoords, |
const int | a_refinementRatio, | ||
const bool | a_interpolateTopography = true , |
||
const bool | a_interpolateThickness = true , |
||
const bool | a_preserveMask = true , |
||
const bool | a_interpolateTopographyGhost = true , |
||
const bool | a_interpolateThicknessGhost = true , |
||
const bool | a_preserveMaskGhost = true , |
||
int | a_thicknessInterpolationMethod = STD_THICKNESS_INTERPOLATION_METHOD |
||
) |
fill the data in this by interpolation from a coarser level
References dx(), exchangeTopography(), getSurfaceHeight(), m_dx, m_floatingMask, m_grids, m_H, m_iceDensity, m_seaLevel, m_topography, m_waterDensity, MAX_THICKNESS_INTERPOLATION_METHOD, IntFineInterp::pwcInterpToFine(), and SMOOTH_SURFACE_THICKNESS_INTERPOLATION_METHOD.
Referenced by MultiLevelDataIBC::initializeIceGeometry(), LevelDataIBC::initializeIceGeometry(), FortranInterfaceIBC::initializeIceGeometry(), MultiLevelDataIBC::regridIceGeometry(), LevelDataIBC::regridIceGeometry(), and FortranInterfaceIBC::regridIceGeometry().
|
inline |
returns true if this object has been defined.
References m_isDefined, and makeCoarser().
LevelSigmaCS * LevelSigmaCS::makeCoarser | ( | int | a_coarsenFactor | ) | const |
returns a new LevelSigmaCS which is a coarsened version of this one (useful for multigrid, etc)
References define(), and LevelSigmaCS().
Referenced by isDefined().
RealVect LevelSigmaCS::mappedCoord | ( | const RealVect & | a_x, |
const DataIndex & | a_index | ||
) | const |
given coordinate in real space, return location in mapped space
given coordinate in real space, return its location in the mapped space
note that a_x is (x,y,z) while mappedCoord returns (sigma,x,y) in 3d. In 2d, mappedCoord returns (x,y)
References m_dx, m_H, and m_topography.
RealVect LevelSigmaCS::realCoord | ( | const RealVect & | a_Xi, |
const DataIndex & | a_index | ||
) | const |
given coordinate in mapped space, return its location in real space note that a_Xi is (sigma,x,y) while realCoord returns (x,y,z) (in 3d – in 2d a_Xi is (x,y) )
given coordinate in mapped space, return its location in real space – this will be a bit slow; probably want to use Fab-based one instead (once it's implemented in Fortran) note that a_Xi is (sigma,x,y) while realCoord returns (x,y,z)
References m_dx, m_H, and m_topography.
Referenced by getNodeRealCoordinates().
void LevelSigmaCS::recomputeGeometry | ( | const LevelSigmaCS * | a_crseCoords, |
const int | a_refRatio | ||
) |
recomputes quantities which are dependent on cell-centered H.
Should be called after cell-centered H is modified to ensure that things don't get out of sync. Also computes face-centered H by averaging to faces
Should be called after H is modified to ensure that things don't get out of sync.
References computeDeltaFactors(), computeSurface(), m_faceH, m_grids, m_H, m_isDefined, and m_topography.
Referenced by IceUtility::eliminateRemoteIce(), getFloatingMask(), AmrIce::incrementIceThickness(), AmrIce::readCheckpointFile(), AmrIce::updateCoordSysWithNewThickness(), and AmrIce::updateGeometry().
void LevelSigmaCS::recomputeGeometryFace | ( | const LevelSigmaCS * | a_crseCoords, |
const int | a_refRatio | ||
) |
recomputes quantities which are dependent on face-centered H.
recomputes quantities which are dependent on cell-centered H.
Should be called after cell-centered H is modified to ensure that things don't get out of sync. Also computes cell-centered H as average of face-centered H
Should be called after H is modified to ensure that things don't get out of sync.
References computeDeltaFactors(), computeSurface(), m_faceH, m_grids, m_H, and m_isDefined.
Referenced by getFloatingMask().
|
inline |
returns sea level
References m_seaLevel.
Referenced by CrevasseCalvingModel::applyCriterion(), DamageConstitutiveRelation::computeFaceMu(), computeSurface(), define(), LevelSigmaCS(), IceUtility::multiplyByGroundedFraction(), and LevelDataIBC::new_thicknessIBC().
|
inline |
set the unit topography difference
for periodic problems : a vector with components (topography(x + Lx , y) - topography(x , y), (topography(x , y + Ly) - topography(x , y) set the background slope for periodic problems entirely on land we need not have periodic topography. Instead we can define a vector slope such that
true topography(x,y) = stored topography(x,y) + slope_x*x + slope_y*y
can be specified. In that case,
surface(x,y) = stored surface(x,y) + slope_x*x + slope_y*y
and
grad(surface) = grad(stored surface(x,y)) + slope
References m_backgroundSlope.
Referenced by defineLevelSigmaCS(), and FortranInterfaceIBC::initializeIceGeometry().
|
protected |
densities of ice and seawater, in kg/(m^3) from Pattyn (2003)
from Vieli and Payne (2005)
References m_backgroundSlope, m_dSigma, m_faceSigma, m_gravity, m_iceDensity, m_seaLevel, m_sigma, and m_waterDensity.
Referenced by LevelSigmaCS().
|
inline |
References m_dSigma, m_faceSigma, and m_sigma.
void LevelSigmaCS::setGradSurface | ( | const LevelData< FArrayBox > & | a_gradSurface | ) |
sets the cell-centered surface gradient
References m_gradSurface.
Referenced by getBackgroundSlope().
void LevelSigmaCS::setGradSurfaceFace | ( | const LevelData< FluxBox > & | a_gradSurface | ) |
sets the face-centred surface gradient
References m_gradSurfaceFace.
Referenced by getBackgroundSlope().
|
inline |
References m_gravity.
|
inline |
interim or not, we also need to set these
References m_iceDensity.
|
inline |
set sea level
References m_seaLevel.
Referenced by HumpIBC::initializeIceGeometry(), VieliPayneIBC::initializeIceGeometry(), and MarineIBC::initializeIceGeometry().
void LevelSigmaCS::setSurfaceHeight | ( | const LevelData< FArrayBox > & | a_surface | ) |
void LevelSigmaCS::setTopography | ( | const LevelData< FArrayBox > & | a_topography | ) |
sets the base height.
In practice, this will probably done by a derived class.
References m_topography.
Referenced by defineLevelSigmaCS(), dx(), HumpIBC::initializeIceGeometry(), VieliPayneIBC::initializeIceGeometry(), MarineIBC::initializeIceGeometry(), and AmrIce::readCheckpointFile().
|
inline |
References m_waterDensity.
void LevelSigmaCS::unitShiftExchange | ( | LevelData< FArrayBox > & | a_level | ) |
References m_Hgrids, and m_unitShift.
|
inline |
seawater density – probably an interim approach
References m_waterDensity.
Referenced by CrevasseCalvingModel::applyCriterion(), DamageConstitutiveRelation::computeFaceMu(), BennCalvingModel::computeRemnant(), VdVCalvingModel::computeRemnant(), computeSurface(), VerticalConductionInternalEnergyIBC::initializeIceInternalEnergy(), L1L2ConstitutiveRelation::modifyTransportCoefficients(), VieliPayneIBC::modifyVelocityRHS(), IceUtility::multiplyByGroundedFraction(), LevelDataIBC::new_thicknessIBC(), VieliPayneBCFunction::operator()(), AmrIce::timeStep(), AmrIce::updateInternalEnergy(), and MarineIBC::velocitySolveBC().
|
protected |
is anything floating?
Referenced by anyFloating(), computeFloatingMask(), and define().
|
protected |
Referenced by computeFloatingMask(), computeSurface(), getBackgroundSlope(), setBackgroundSlope(), and setDefaultValues().
|
protected |
Referenced by computeDeltaFactors(), define(), deltaFactors(), and LevelSigmaCS().
|
protected |
Referenced by define(), getDSigma(), LevelSigmaCS(), setDefaultValues(), and setFaceSigma().
|
protected |
Referenced by computeDeltaFactors(), computeSurface(), define(), dx(), getNodeRealCoordinates(), interpFromCoarse(), mappedCoord(), and realCoord().
|
protected |
Referenced by define(), faceDeltaFactors(), and LevelSigmaCS().
|
protected |
face-centered ice thickness
Referenced by define(), getFaceH(), LevelSigmaCS(), recomputeGeometry(), and recomputeGeometryFace().
|
protected |
Referenced by define(), getFaceSigma(), LevelSigmaCS(), setDefaultValues(), and setFaceSigma().
|
protected |
floating mask
Referenced by computeFloatingMask(), computeSurface(), define(), getFloatingMask(), and interpFromCoarse().
|
protected |
ghost vector
Referenced by define(), and ghostVect().
|
protected |
cell-centered gradient of surface elevation
Referenced by computeSurface(), define(), getGradSurface(), LevelSigmaCS(), and setGradSurface().
|
protected |
face-centered gradient of surface elevation
Referenced by computeSurface(), define(), getGradSurfaceFace(), LevelSigmaCS(), and setGradSurfaceFace().
|
protected |
acceleration due to gravity – probably an interim approach
Referenced by gravity(), setDefaultValues(), and setGravity().
|
protected |
|
protected |
cell-centered ice thickness
Referenced by computeDeltaFactors(), computeFloatingMask(), computeSurface(), define(), getH(), getSurfaceHeight(), interpFromCoarse(), LevelSigmaCS(), mappedCoord(), realCoord(), recomputeGeometry(), and recomputeGeometryFace().
|
protected |
"flattened" grids (in 3d) for quantities like H
Referenced by define(), and unitShiftExchange().
|
protected |
ice density
Referenced by computeFloatingMask(), getSurfaceHeight(), iceDensity(), interpFromCoarse(), setDefaultValues(), and setIceDensity().
|
protected |
has this object been defined?
Referenced by define(), isDefined(), LevelSigmaCS(), recomputeGeometry(), and recomputeGeometryFace().
|
protected |
Referenced by computeFloatingMask(), computeSurface(), define(), getSurfaceHeight(), interpFromCoarse(), LevelSigmaCS(), seaLevel(), setDefaultValues(), and setSeaLevel().
|
protected |
Referenced by define(), getSigma(), LevelSigmaCS(), setDefaultValues(), and setFaceSigma().
|
protected |
cell-centered surface elevation
Referenced by computeSurface(), define(), getSurfaceHeight(), LevelSigmaCS(), and setSurfaceHeight().
|
protected |
cell-centered thickness over flotation
Referenced by computeSurface(), define(), getThicknessOverFlotation(), and LevelSigmaCS().
|
protected |
cell-centered topography
Referenced by computeDeltaFactors(), computeFloatingMask(), computeSurface(), define(), exchangeTopography(), getSurfaceHeight(), getTopography(), interpFromCoarse(), LevelSigmaCS(), mappedCoord(), realCoord(), recomputeGeometry(), and setTopography().
|
protected |
Referenced by unitShiftExchange().
|
protected |
seawater density – probably an interim approach
Referenced by computeFloatingMask(), getSurfaceHeight(), interpFromCoarse(), setDefaultValues(), setWaterDensity(), and waterDensity().