Mushy Layer  1.0
Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
AMRLevelMushyLayer Class Reference

AMRLevel for mushy layer calculations. More...

#include <AMRLevelMushyLayer.H>

Inheritance diagram for AMRLevelMushyLayer:
AMRLevel

Public Member Functions

 AMRLevelMushyLayer ()
 Default constructor.
 
 AMRLevelMushyLayer (MushyLayerOptions a_opt, MushyLayerParams a_params)
 Full constructor. Arguments are same as in define()
 
void define (MushyLayerOptions a_opt, MushyLayerParams a_params)
 Defines this AMRLevelMushyLayer. More...
 
virtual ~AMRLevelMushyLayer ()
 Destructor.
 
virtual void define (AMRLevel *a_coarserLevelPtr, const Box &a_problemDomain, int a_level, int a_refRatio)
 Never called: historical.
 
void setDefaults ()
 Set default values. More...
 
virtual void define (AMRLevel *a_coarserLevelPtr, const ProblemDomain &a_problemDomain, int a_level, int a_refRatio)
 Define new AMRLevelMushyLayer from coarser. More...
 
virtual Real advance ()
 Advance by one timestep.
 
virtual void postTimeStep ()
 Things to do after a timestep.
 
void doExplicitReflux (int a_var)
 Explicitly reflux scalar field.
 
Real doHCreflux ()
 Reflux enthalpy and bulk concentration.
 
void doMomentumReflux (Vector< LevelData< FArrayBox > *> &compVel)
 Do reflux on momentum field.
 
virtual void tagCells (IntVectSet &a_tags)
 Create tags for regridding.
 
void tagMushLiquidBoundary (IntVectSet &localTags)
 Tags cells on the mush side of the mush-liquid boundary.
 
void tagBoundaryLayerCells (IntVectSet &localTags)
 Tag cells on the domain boundary.
 
void tagCenterCells (IntVectSet &localTags, Real radius=0)
 Tag cells inside grid.
 
void tagMushyCells (IntVectSet &localTags)
 Tag cells with porosity < 1.
 
void tagCellsVar (IntVectSet &localTags, Real refineThresh, int taggingVar, int vectorTaggingVar, TaggingMethod taggingMethod, int comp=-1)
 Add cells that specify criteria.
 
void levelSetup ()
 Setup menagerie of data structures.
 
void defineCFInterp ()
 Setup CF interp objects.
 
void setEta (Real a_eta)
 Set eta on all projection operators. More...
 
Real maxAllowedEta ()
 Maximum allowed eta.
 
void setAdvVelCentering (Real a_fraction)
 Set adv_vel_centering on all levels.
 
virtual void tagCellsInit (IntVectSet &a_tags)
 Create tags at initialization.
 
virtual void regrid (const Vector< Box > &a_newGrids)
 Set up data on this level after regridding.
 
virtual void initialGrid (const Vector< Box > &a_newGrids)
 Initialize grids.
 
virtual void postInitialGrid (const bool a_restart)
 Operations to do after initial grids have been set up.
 
virtual void postRegrid (int a_base_level)
 Operations to execute after we've done regridding.
 
void createDataStructures ()
 Create structures to hold data.
 
void initDataStructures ()
 Initialize data structures - resize arrays etc.
 
void computeAllVelocities (bool doFRupdates)
 Compute cell-centred and face-centred velocities.
 
virtual void initialData ()
 Initialize data.
 
void initialDataCornerFlow ()
 Initial data for corner flow problem.
 
void initialDataHRL ()
 Initial data for Horton-Rogers-Lapwood problem.
 
void initialDataConvectionMixedPorous ()
 Initial data for convection in a hetrogeneous porous medium.
 
void initialDataPoiseuille ()
 Initial data for poiseuille flow.
 
void initialDataIceBlock ()
 Initial data for a solid ice block.
 
void initialDataMushyLayer ()
 Initial data for mushy layer simulations.
 
void initialDataPorousHole ()
 Mushy layer with a porous hole in the middle.
 
void initialDataZeroPorosityTest ()
 Initial data for a test case with regions of zero porosity.
 
void initialDataDiffusion ()
 Initial data for diffusion problem.
 
void initialDataVortexPair ()
 Initial data for vortex pair problem.
 
void addVortex (RealVect center, Real strength, Real radius)
 Add vortex to domain.
 
void initialDataRayleighBenard ()
 Initial data for Rayleigh-Benard problem.
 
void initialDataSoluteFlux ()
 Initial data for solute flux test problem.
 
void initialDataSidewallHeating ()
 Initial data for sidewall heating problem.
 
void initialDataDefault ()
 Default initial data.
 
void initialDataBurgers ()
 Initial data for burgers problem.
 
void initialDataRefluxTest ()
 Initial data for problem to test refluxing.
 
bool loadAdvVel ()
 Whether or not we're loading advection velocities.
 
virtual bool convergedToSteadyState ()
 Have we converged to steady state?
 
virtual void postInitialize ()
 Things to do after initialization.
 
void computeVorticityStreamfunction ()
 Calculate vorticity and streamfunction so we can plot them.
 
void computeVorticity ()
 Compute vorticity.
 
void getExtraPlotFields ()
 Get $ \pi $ and $ \phi $ from projection so we can plot them.
 
void writeAMRHierarchy (string filename)
 Write AMR hierarchy to file.
 
void setDimensionlessReferenceInitial ()
 Set reference point for nondimensionalisation to be the initial temperature/salinity.
 
void setDimensionlessReferenceEutectic ()
 Set reference point for nondimensionalisation to be the eutectic temperature/salinity.
 
virtual Real computeDt ()
 Returns the dt computed earlier for this level.
 
virtual Real computeInitialDt ()
 Compute dt using initial data.
 
void computeDiagnostics ()
 Compute nusselt number, salt fluxes etc.
 
void defineSolvers (Real a_time)
 Define solvers.
 
void computeChimneyDiagnostics ()
 Compute chimney diagnostics.
 
DisjointBoxLayout grids ()
 Whether or not to compute diagnostics. More...
 
void shiftData (int dir, int distance)
 Shift data in domain.
 
void dx (Real newDx)
 Set $ \Delta x $ on this level.
 
Real dx ()
 Get $ \Delta x $ on this level.
 
void reshapeData (DisjointBoxLayout newGrids, ProblemDomain newDomain)
 Transfer data to new disjoint box layout. More...
 
void refine (Real ref_ratio, DisjointBoxLayout a_grids, ProblemDomain a_domain)
 Refine data.
 
void doPostRegridSmoothing (bool a_smoothVel=true, bool a_smoothScalar=true)
 Do smoothing after regridding.
 
void setSmoothingCoeff (Real a_coeff)
 Set the coefficient for smoothing fields.
 
void setSmoothingDone (bool a_smoothingDone)
 
void addMeltPond (int depth, Real salinity, Real enthalpy, bool rescaleExistingSolution=false)
 Add a melt pond at the top of the domain. More...
 
void setVelZero (LevelData< FArrayBox > &a_vel, Real a_limit=-1, int a_radius=0)
 Make sure zero porosity regions have no velocity (cell centred velocity)
 
void setVelZero (LevelData< FluxBox > &a_vel, Real a_limit=-1)
 Make sure zero porosity regions have no velocity (face centred velocity)
 
void setCCVelZero (Real a_limit)
 Set velocity to zero when porosity is below specified limit (cell centred)
 
void updateEnthalpyVariables ()
 Update thermodynamic variables based on phase diagram.
 
void computeLambdaPorosity ()
 Compute lambda/porosity.
 
void set_compute_diagnostics (bool compute_diags)
 Set whether or not to compute diagnostics.
 
void horizAverage ()
 Average enthalpy and bulk salinity then recalculate all other fields.
 
void smoothScalarField (LevelData< FArrayBox > &a_phi, int a_var, Real a_smoothing)
 Smooth some scalar field. More...
 
void smoothEnthalpyBulkConc (Real a_smoothing)
 
- Public Member Functions inherited from AMRLevel
virtual void preRegrid (int a_base_level, const Vector< Vector< Box > > &a_new_grids)
 
virtual void conclude (int a_step) const
 
virtual void writeCheckpointHeader (HDF5Handle &a_handle) const=0
 
virtual void writeCheckpointLevel (HDF5Handle &a_handle) const=0
 
virtual void readCheckpointHeader (HDF5Handle &a_handle)=0
 
virtual void readCheckpointLevel (HDF5Handle &a_handle)=0
 
virtual void writePlotHeader (HDF5Handle &a_handle) const=0
 
virtual void writePlotLevel (HDF5Handle &a_handle) const=0
 
virtual void writeCustomPlotFile (const std::string &a_prefix, int a_step) const
 
virtual void finerLevelPtr (AMRLevel *a_finer_level_ptr)
 
virtual void dt (Real a_dt)
 
virtual void time (Real a_time)
 
virtual void initialDtMultiplier (Real a_initial_dt_multiplier)
 
virtual Real dt () const
 
virtual Real time () const
 
virtual Real initialDtMultiplier () const
 
virtual const ProblemDomainproblemDomain () const
 
virtual Vector< Boxboxes () const
 
bool isDefined () const
 
bool hasCoarserLevel () const
 
bool hasFinerLevel () const
 
virtual int level () const
 
virtual int refRatio () const
 
Vector< AMRLevel *> getAMRLevelHierarchy ()
 

Public Attributes

Diagnostics m_diagnostics
 Diagnostics e.g. Nusselt number, solute flux.
 

Protected Types

enum  solverRestart { m_restartHalveDt, m_restartResetData }
 What to do when restarting after a failed timestep. More...
 

Protected Member Functions

Real convergedToSteadyState (const int a_var, bool vector=false)
 Has a particular variable converged to steady state?
 
void compute_d_dt (const int a_var, LevelData< FArrayBox > &diff, bool vector=false)
 Compute $ \partial / \partial t $ of some field.
 
bool solvingFullDarcyBrinkman ()
 Determine if we're solving the full time dependent equation, or just darcy's law.
 
void updateEnthalpyVariablesOld ()
 Update thermodynamic variables based on phase diagram (old time)
 
bool finestLevel ()
 Returns true if there is not a finer level with a grid defined on it.
 
int getFinestLevel ()
 Return the finest level.
 
Real computeDt (Real cfl)
 Returns $ \Delta t $.
 
Real getMaxVelocity ()
 Get largest velocity in domain.
 
Real getMaxVelocityForCFL ()
 get largest advection velocity, considering $ \mathbf{U}/\chi $ if appropriate
 
Real computeMaxUChi ()
 Compute max $ \mathbf{U}/\chi $.
 
Real getMaxAdvVel ()
 Get max advection velocity.
 
Real computeDt (bool growdt)
 Returns $ \Delta t $.
 
void restartTimestepFromBackup (bool ignorePressure=false)
 Replace with data from before the last backup.
 
void backupTimestep ()
 Save the current data.
 
void getHierarchyAndGrids (Vector< AMRLevelMushyLayer *> &a_hierarchy, Vector< DisjointBoxLayout > &a_grids, Vector< int > &a_refRat, ProblemDomain &a_lev0Dom, Real &a_lev0Dx)
 Get information on the entire AMR hierarchy.
 
void calculatePermeability ()
 Calculate permeability $ \Pi $.
 
void calculateCoarseFineBoundaries (int a_var, bool vector=false)
 Fills coarse-fine boundaries on variable by quadratic interpolation.
 
int multiCompAdvectDiffuse (LevelData< FArrayBox > &a_phi_old, LevelData< FArrayBox > &a_phi_new, LevelData< FArrayBox > &a_src, bool doFRupdates=false, bool computeAdvectiveSrc=true)
 Advection diffusion for multiple components.
 
void computeScalarAdvectiveSrc (LevelData< FArrayBox > &a_src, int a_var, bool converged, int a_comp=0)
 Compute src term and add it to a_src.
 
void computeScalarAdvectiveSrcHC (LevelData< FArrayBox > &a_src, LevelData< FluxBox > &edgeScalTotal, bool converged)
 Compute multi component source term: $ (\mathbf{U} \cdot \nabla T, \mathbf{U} \cdot \nabla S_l)$.
 
void incrementHCFluxRegisters (LevelData< FluxBox > &flux, Real fluxMult)
 Increment Enthalpy-Bulk concentration flux registers with flux.
 
void computeTotalAdvectiveFluxes (LevelData< FluxBox > &edgeScalTotal)
 Compute total fluxes that provide the source term for implicit enthalpy-salinity updates. More...
 
void advectScalar (const int a_scalarVar, const int a_advectionVar, LevelData< FluxBox > &a_advVel, bool doFRupdates=true)
 Update a scalar field by advection. More...
 
void advectScalar (const int a_scalarVar, const int a_advectionVar, LevelData< FluxBox > &a_advVel, bool doFRupdates, LevelData< FluxBox > &flux)
 Update a scalar field with advective terms.
 
void updateScalarFluxRegister (int a_scalarVar, LevelData< FluxBox > &flux, Real scale)
 Update flux register for scalar field.
 
void advectLambda (bool doFRupdates=true)
 Advect lambda field (for freestream preservation)
 
void computeAdvectionVelSourceTerm (LevelData< FArrayBox > &advectionSourceTerm)
 Make source term for the advection velocity solve.
 
void horizontallyAverage (LevelData< FArrayBox > &a_averaged, LevelData< FluxBox > &a_phi)
 Average vertical fluxes in the horizontal direction.
 
void horizontallyAverage (Vector< Real > &averageVector, LevelData< FluxBox > &a_phi)
 Average vertical fluxes in the horizontal direction.
 
void horizontallyAverage (LevelData< FArrayBox > &a_averaged, LevelData< FluxBox > &a_phi, Vector< Real > &averageVector)
 Average vertical fluxes in the horizontal direction.
 
void horizontallyAverage (LevelData< FArrayBox > &a_averaged, LevelData< FArrayBox > &a_phi)
 Average cell centred fields in the horizontal direction.
 
void horizontallyAverage (Vector< Real > &averageVector, LevelData< FArrayBox > &a_phi)
 Average cell centred fields in the horizontal direction.
 
void horizontallyAverage (LevelData< FArrayBox > &a_averaged, LevelData< FArrayBox > &a_phi, Vector< Real > &globalAveraged)
 Average cell centred values in the horizontal direction.
 
Real averageOverLiquidRegion (int a_var)
 Compute average value of some scalar field over the liquid region of the domain.
 
void computeRadianceIntensity ()
 Compute the intensity of light from the surface at each point in the domain.
 
void advectPassiveTracer ()
 Advect a passive tracer $ \xi $. More...
 
void computeScalarDiffusiveSrc (int a_scalarBulkConc, LevelData< FArrayBox > &a_src)
 Compute diffusive src for some scalar with bulk concentration $ \xi $. More...
 
void computeScalarConcInLiquid (LevelData< FArrayBox > &liquid_tracer_conc, int a_tracerVar)
 Compute concentration of scalar (with bulk concentration $ \xi $ ) in the liquid phase. More...
 
void advectTracer (int a_tracerVar, LevelData< FArrayBox > &a_src)
 Advect generic tracer.
 
void advectActiveTracer ()
 Advect active tracer $ \zeta $. More...
 
void computeActiveTracerSourceSink (LevelData< FArrayBox > &a_srcSink)
 Compute sources/sinks for active tracer.
 
Real computeMushDepth (Real a_porosity_criteria=0.99)
 Compute mushy layer depth. More...
 
void getTotalFlux (LevelData< FluxBox > &totalFlux)
 Calculate advective + diffusive flux.
 
void computeLapVel (LevelData< FArrayBox > &a_lapVel, LevelData< FArrayBox > &a_vel, const LevelData< FArrayBox > *a_crseVelPtr)
 Compute $ \nabla^2 \mathbf{U} $.
 
void computeScalDiffusion (const int a_var, LevelData< FArrayBox > &a_lapScal, Real a_time)
 Compute $ \nabla^2 \phi $.
 
void computeScalDiffusion (LevelData< FArrayBox > &diffusiveSrc, Real a_time)
 Compute $ \nabla^2 (H, S) $.
 
void computeScalDiffusion (LevelData< FArrayBox > &a_src, int a_var)
 Compute $ \nabla^2 \psi $ for some scalar field $ \psi $.
 
void computeCCvelocity (const LevelData< FArrayBox > &advectionSourceTerm, Real a_oldTime, Real a_dt, bool doFRupdates=false, bool a_doProjection=true, bool compute_uDelU=true, bool a_MACprojection=false)
 Compute new cell centred velocity. More...
 
bool isVelocityTimeDependent ()
 Is velocity time dependent?
 
void computeUstar (LevelData< FArrayBox > &a_uDelu, const LevelData< FArrayBox > &advectionSourceTerm, Real a_oldTime, Real a_dt, bool doFRupdates, bool a_MACprojection=false, bool compute_uDelU=true)
 Compute $ \mathbf{U}^* $, the velocity before we project it. More...
 
void computeUstarSrc (LevelData< FArrayBox > &src, const LevelData< FArrayBox > &advectionSourceTerm, Real src_time, bool a_MACprojection=false, bool compute_uDelU=true)
 Compute the source term for the $\mathbf{U}^* $ solve.
 
void computeUDelU (LevelData< FArrayBox > &U_adv_src, const LevelData< FArrayBox > &advectionSourceTerm, Real a_oldTime, Real a_dt)
 Compute $ \mathbf{u} \cdot \nabla \left( \mathbf{U}/\chi \right) $.
 
bool isCurrentCFLSafe (bool printWarning=false)
 Check CFL based on latest velocity is safe to use for explicit updates.
 
void addHeatSource (LevelData< FArrayBox > &src)
 Add a fixed heat source. More...
 
void setVelZero (FArrayBox &a_vel, const FArrayBox &a_porosity, const Real a_limit, const int a_radius=0)
 Set velocity to zero where porosity is sufficiently small (FArrayBox version)
 
void computeGradP (LevelData< FArrayBox > &gradP, Real a_time, bool a_macProjection=false)
 Returns $ \nabla P $ at the specified time.
 
void traceAdvectionVel (LevelData< FluxBox > &a_advVel, LevelData< FArrayBox > &a_old_vel, LevelData< FArrayBox > &U_chi, const LevelData< FArrayBox > &a_viscousSource, PatchGodunov &a_patchGodVelocity, Real a_old_time, Real a_dt)
 Trace the advection velocity forward in time.
 
void computePredictedVelocities (LevelData< FluxBox > &U_chi_new, LevelData< FArrayBox > &a_traceVel, LevelData< FluxBox > &a_advVel, LevelData< FArrayBox > &U_chi, const LevelData< FArrayBox > &a_viscousSource, PatchGodunov &a_patchGodVelocity, LevelData< FluxBox > &a_grad_eLambda, LevelData< FluxBox > &a_gradPhi, LevelData< FluxBox > &porosityFace, Real a_old_time, Real a_dt)
 Compute advection velocities at the half time step.
 
void computeScalarAdvectiveFlux (LevelData< FluxBox > &a_edgeScal, LevelData< FArrayBox > &a_old_scal, LevelData< FluxBox > &a_adv_vel, LevelData< FluxBox > &a_inflowOutflowVel, LevelData< FArrayBox > &a_old_vel, LevelData< FArrayBox > &a_diffusiveSrc, PatchGodunov &a_patchGod, Real a_old_time, Real a_dt)
 Compute the advective flux given the scalar field, advection velocity and source term specified. More...
 
void upwind (LevelData< FluxBox > &a_edgeScal, LevelData< FArrayBox > &a_old_scal, LevelData< FluxBox > &a_adv_vel, LevelData< FluxBox > &a_inflowOutflowVel, LevelData< FArrayBox > &a_old_vel, LevelData< FArrayBox > &a_diffusiveSrc, PatchGodunov &a_patchGod, Real a_old_time, Real a_dt)
 Refactored this so we can also use for velocities if we want.
 
void computeScalarAdvectiveFlux (LevelData< FluxBox > &a_edgeScal, LevelData< FArrayBox > &a_scalar_advection_old, LevelData< FArrayBox > &a_src, LevelData< FluxBox > &a_advVel, int a_advectionVar, Real a_old_time, Real a_dt)
 Compute advective flux.
 
void computeScalarAdvectiveFlux (LevelData< FluxBox > &a_edgeScal, int a_advectionVar, int a_diffusionVar, LevelData< FluxBox > &a_advVel, Real a_old_time, Real a_dt)
 Compute the advective flux. More...
 
void computeScalarAdvectiveFluxMultiComp (LevelData< FluxBox > &a_edgeScal, LevelData< FluxBox > &a_advVel, PatchGodunov &a_patchGod, LevelData< FArrayBox > &a_scalOld, Real a_old_time, Real a_dt)
 Compute the advective flux for the enthalpy and salinity equations,. More...
 
void predictVelocities (LevelData< FArrayBox > &a_uDelU, LevelData< FluxBox > &a_advVel, const LevelData< FArrayBox > &a_src, Real old_time, Real a_dt, bool doFRupdates)
 Calculate $ \mathbf{u} \cdot \nabla \left( \mathbf{u}/\chi \right) $.
 
void calculateTimeIndAdvectionVel (Real a_time, LevelData< FluxBox > &a_advVel)
 Calculate advection velocity for momentum equations which are time independent. More...
 
void fillUnprojectedDarcyVelocity (LevelData< FluxBox > &a_advVel, Real time)
 Fill unprojected darcy velocity. More...
 
int getMaxLevel ()
 Get maximum level allowed in AMR hierarchy.
 
void computeViscosity ()
 Compute the viscosity given the current liquid concentration.
 
void fillAdvVel (Real time, LevelData< FluxBox > &a_advVel)
 Fill advection velocity ghost cells.
 
void fillAnalyticVel (LevelData< FluxBox > &a_advVel)
 Fill the velocity with some analytically determined function.
 
void fillAnalyticVel (FArrayBox &velDir, int dir, int comp, bool project)
 Fill the velocity with some analytically determined function.
 
void fillAnalyticVel (LevelData< FArrayBox > &velDir)
 Fill the velocity with some analytically determined function.
 
void fillBuoyancy (LevelData< FArrayBox > &a_buoyancy, Real a_time)
 Compute the buoyancy force.
 
void fillBuoyancy (LevelData< FArrayBox > &a_buoyancy, LevelData< FArrayBox > &a_temperature, LevelData< FArrayBox > &a_liquidConc, LevelData< FArrayBox > &a_porosity)
 Compute the buoyancy force.
 
void fillBuoyancy (FArrayBox &buoyancy, FArrayBox &temperature, FArrayBox &liquidConc, FArrayBox &porosity, FArrayBox &bodyForce)
 Compute the buoyancy force. More...
 
void fillPressureSrcTerm (LevelData< FArrayBox > &gradP, LevelData< FArrayBox > &pressureScale, Real a_time, bool a_MACprojection)
 Compute $ \nabla P $ at time $ t $.
 
void fillAMRVelPorosity (Vector< LevelData< FArrayBox > *> &amrVel, Vector< RefCountedPtr< LevelData< FluxBox > > > &amrPorosityFace, Vector< RefCountedPtr< LevelData< FArrayBox > > > &amrPorosity)
 Compile velocity and porosity across an AMR hierarchy.
 
void fillHC (LevelData< FArrayBox > &a_phi, Real a_time, bool doInterior=true, bool quadInterp=true)
 Should include correct exchange calls.
 
void fillFixedPorosity (LevelData< FArrayBox > &a_porosity)
 Fill porosity with some fixed value if appropriate.
 
void fillTCl (LevelData< FArrayBox > &a_phi, Real a_time, bool doInterior=true, bool quadInterp=true)
 Should include correct exchange calls.
 
void fillMultiComp (LevelData< FArrayBox > &a_phi, Real a_time, int scal1, int scal2, bool doInterior=true, bool quadInterp=true, bool apply_bcs=true)
 Should include correct exchange calls.
 
void fillAMRLambda (Vector< LevelData< FArrayBox > *> &amrLambda)
 Compile $ \lambda $ across the entire AMR hierarchy.
 
void calculateAnalyticSolns (bool enforceSolutions=true)
 Calculate analytic solution to this problem (if possible)
 
void setPorosity (FArrayBox &a_porosity)
 Set fixed porosity.
 
void fillFrameVelocity ()
 Fill frame advection velocity.
 
void addPerturbation (int a_var, Real alpha, int waveNumber=-1, Real phaseShift=0)
 Add a small perturbation to the solution.
 
void addMeltPond ()
 Add a melt pond if we want.
 
void defineUstarSolver (Vector< RefCountedPtr< LevelBackwardEuler > > &a_UstarBE, Vector< RefCountedPtr< LevelTGA > > &a_UstarTGA)
 Define implicit solver for $ \mathbf{u}^* $ including timestepping.
 
void defineUstarMultigrid ()
 Define multigrid solver for $ \mathbf{u}^* $.
 
void initializeLevelPressure (Real a_currentTime, Real a_dtInit)
 Initialize pressure on this level.
 
void initializeGlobalPressure (Real dtInit=-1, bool init=true)
 Initialize pressure across all levels.
 
void computeInitAdvectionVel ()
 Compute initial advection velocities.
 
void initTimeIndependentPressure (AMRLevelMushyLayer *lev, int a_max_num_iter=-1)
 Initialise phi (pressure for Darcy equation) More...
 
void resetLambda ()
 Flag if a new level has been added. More...
 
void AMRRefluxLambda ()
 Reflux $ \lambda $ across AMR hierarchy. More...
 
void defineRegridAMROp (AMRPoissonOpFactory &a_factory, const Vector< DisjointBoxLayout > &a_grids, const Vector< ProblemDomain > &a_domains, const Vector< Real > &a_amrDx, const Vector< int > &a_refRatios, const int &a_lBase)
 Define operator for post regrid smoothing.
 
void getCoarseScalarDataPointers (const int a_scalarVar, LevelData< FArrayBox > **a_coarserDataOldPtr, LevelData< FArrayBox > **a_coarserDataNewPtr, LevelFluxRegister **a_coarserFRPtr, LevelFluxRegister **a_finerFRPtr, Real &a_tCoarserOld, Real &a_tCoarserNew)
 Smooth the velocity field. More...
 
void fillVectorField (LevelData< FArrayBox > &a_vector, Real a_time, int a_var, bool doInterior=false, bool quadInterp=false)
 Doesn't do any time interpolation - just fills ghost cells if necessary.
 
void fillScalars (LevelData< FArrayBox > &a_scal, Real a_time, const int a_var, bool doInterior=false, bool quadInterp=false, int a_comp=0, bool apply_bcs=true)
 Fill ghost cells of scalar field.
 
void fillBuoyancy (LevelData< FArrayBox > &a_buoyancy, Real a_time, bool quadInterp)
 Helper function.
 
void fillScalarFace (LevelData< FluxBox > &a_scal, Real a_time, const int a_var, bool doInterior=false, bool quadInterp=false)
 Fill ghost cells of scalar field (edge-centred)
 
void fillScalarFace (LevelData< FluxBox > &a_scal, Real a_time, const int a_var, CellToEdgeAveragingMethod method, bool doInterior=false, bool quadInterp=false, Real smoothing=0.0)
 Fill ghost cells of scalar field (edge-centred)
 
void getScalarBCs (BCHolder &thisBC, int a_var, bool a_homogeneous)
 Get boundary conditions for implicit problems.
 
void computeInflowOutflowAdvVel ()
 Refactored this to be consistent.
 
PhysIBCgetScalarIBCs (int a_var)
 Get initial and boundary conditions for advection.
 
void convertBCType (const int a_implicitBC, const Real a_implicitVal, int &a_explicitBC, Real a_explicitVal)
 Need to convert implicit dirichlet/neuman/inflow outflow BCs to advection BCs.
 
int convertBCType (const int a_implicitBC)
 Convert implicit BCs.
 
bool crashed ()
 Check the simulation hasn't crashed.
 
void doRegularisationOps (LevelData< FArrayBox > &a_scal, int a_var, int a_comp=0)
 Operations to keep variables bounded correctly. More...
 
void doRegularisationOps (LevelData< FluxBox > &a_scal, int a_var, int a_comp=0)
 Operations to keep variables bounded correctly.
 
void doRegularisationOps (int a_var, FArrayBox &a_state, int a_comp=0)
 Operations to keep variables bounded correctly.
 
void doRegularisationOpsNew (int a_var, FArrayBox &a_state, int a_comp=0)
 New version of doRegularisationOps() which uses fortran routines.
 
void setFluxRegistersZero ()
 Reset all flux registers.
 
void copyNewToOldStates ()
 Copy new data to old data holders.
 
bool doVelocityAdvection ()
 Do we need to do velocity advection? More...
 
void computeAdvectionVelocities (LevelData< FArrayBox > &advectionSrc, Real advVelCentering=0.5)
 Compute advection velocity.
 
void correctEdgeCentredVelocity (LevelData< FluxBox > &a_advVel, Real a_dt)
 Correct velocity (projection + volume discrepancy)
 
Real computeDtInit (int finest_level)
 Compute initial $ \Delta t $.
 
AMRLevelMushyLayergetCoarserLevel () const
 Get the next coarser level.
 
AMRLevelMushyLayergetCoarsestLevel ()
 Get the coarsest level.
 
AMRLevelMushyLayergetFinerLevel () const
 Get the next finer level.
 
void setAdvectionBCs ()
 Set BCs on advection physics classes.
 
void defineIBCs ()
 Define initial and boundary conditions for advection problems.
 

Protected Attributes

LevelData< FluxBoxm_advVel
 Advection velocity.
 
LevelData< FluxBoxm_advVelOld
 Advection velocity at old time.
 
LevelData< FluxBoxm_advVelNew
 Advection velocity at new time.
 
LevelData< FluxBoxm_frameAdvVel
 Frame advection velocity.
 
LevelData< FluxBoxm_totalAdvVel
 Frame advection velocity + fluid advection velocity.
 
Real AMRSaltSum_new
 Total salt content in domain at new time.
 
Real AMRSaltSum_old
 Total salt content in domain at old time.
 
Real AMREnthalpySum_new
 Total heat content in domain at new time.
 
Real AMREnthalpySum_old
 Total heat content in domain at new time.
 
LevelData< FArrayBoxm_saltFluxTop
 Salt flux at the top of the domain.
 
LevelData< FArrayBoxm_saltFluxBottom
 Salt flux at the bottom of the domain.
 
LevelData< FArrayBoxm_frameVel
 Frame velocity (cell centred). Don't think we need this.
 
Vector< int > m_outputScalarVars
 List of scalar vars to write out to HDF5 files.
 
Vector< int > m_outputVectorVars
 List of vector vars to write out to HDF5 files.
 
Vector< int > m_chkVectorVars
 List of vector vars to write to checkpoint files.
 
Vector< int > m_chkScalarVars
 List of scalar vars to write to checkpoint files.
 
int m_numOutputComps
 Number of component we write out to HDF5 files.
 
Real m_maxLambda
 Store max value of $ \lambda-1 $ from previous timestep. More...
 
Vector< RefCountedPtr< LevelData< FArrayBox > > > m_scalarNew
 Scalar fields at new time.
 
Vector< RefCountedPtr< LevelData< FArrayBox > > > m_scalarOld
 Scalar fields at old time.
 
Vector< RefCountedPtr< LevelData< FArrayBox > > > m_scalarRestart
 Backed up scalar fields to be used for restarting.
 
LevelData< FArrayBoxm_dPorosity_dt
 Lagged $ \partial \chi / \partial t $.
 
Vector< RefCountedPtr< LevelData< FArrayBox > > > m_vectorNew
 Vector fields at new time.
 
Vector< RefCountedPtr< LevelData< FArrayBox > > > m_vectorOld
 Vector fields at old time.
 
Vector< RefCountedPtr< LevelData< FArrayBox > > > m_vectorRestart
 Backed up vector fields to be used for restarting.
 
Vector< string > m_scalarVarNames
 Names of scalar variables.
 
Vector< string > m_vectorVarNames
 Names of vector variables.
 
Vector< int > m_vectRestartVars
 Indexes of vector fields to save for restarting.
 
Vector< int > m_scalRestartVars
 Indexes of scalar fields to save for restarting.
 
RefCountedPtr< LevelTGAm_enthalpySalinityTGA
 TGA integrators (2nd order in time)
 
RefCountedPtr< LevelBackwardEulerm_enthalpySalinityBE
 Backward Euler integrators (1st order in time)
 
RefCountedPtr< AMRMultiGrid< LevelData< FArrayBox > > > m_diffuseAMRMG [m_numScalarVars]
 AMRMultigrid operators for each scalar.
 
RefCountedPtr< AMRFASMultiGrid< LevelData< FArrayBox > > > m_diffuseAMRFASMG [m_numScalarVars]
 FASAMRMultigrid operators for each scalar.
 
RefCountedPtr< AMRFASMultiGrid< LevelData< FArrayBox > > > m_multiCompFASMG
 Multicomponent FAS Multigrid solver (for enthalpy-bulk concentration solves)
 
RefCountedPtr< AMRLevelOpFactory< LevelData< FArrayBox > > > m_HCOpFact
 Operator factories for scalar vars.
 
RefCountedPtr< AMRMultiGrid< LevelData< FArrayBox > > > m_uStarAMRMG [SpaceDim]
 AMRMultigrid solver for unprojected velocity, $ \mathbf{u}^* $.
 
RefCountedPtr< AMRLevelOpFactory< LevelData< FArrayBox > > > m_uStarOpFact [SpaceDim]
 Operator factory for unprojected velocity, $ \mathbf{u}^* $.
 
RefCountedPtr< AMRPoissonOpm_viscousOp [SpaceDim]
 Operator for computing components of $ \nabla^2 \mathbf{u} $.
 
bool m_isDefined
 Is object defined.
 
MushyLayerOptions m_opt
 Struct to contain all the options for this simulation. More...
 
Real m_computedCFL
 The CFL number we're currently running at. More...
 
Real m_domainWidth
 Domain width.
 
Real m_domainHeight
 Domain height.
 
IntVect m_numCells
 Number of cells in each direction of domain.
 
Real m_max_dt_growth
 Maxmimum fractional growth in $ \Delta t $ between timesteps.
 
bool m_useLimiting
 Should we do slope limiting in advection solves?
 
Real m_dx
 Cell spacing.
 
Real m_scalarDiffusionCoeffs [m_numScalarVars]
 Diffusion coefficients for different scalar fields.
 
MushyLayerParams m_parameters
 Physical parameters.
 
Real m_prev_diag_output
 Time when we last produced diagnostics output.
 
bool m_enforceGradP
 Should we enforce some analytic $ \nabla P $. More...
 
bool m_usePrevPressureForUStar
 Should we use previous $ \nabla P $ to calculate $ \mathbf{u}^* $. More...
 
int m_pressureScaleVar
 Which field scales the pressure. More...
 
PhysBCUtilm_physBCPtr
 Pointer to physics boundary condition object.
 
FineInterp m_fineInterpScalar
 For interpolating scalars from coarser grids.
 
FineInterp m_fineInterpVector
 For interpolating vectors from coarser grids.
 
CoarseAverage m_coarseAverageScalar
 For averaging scalars to coarse grids.
 
CoarseAverage m_coarseAverageHC
 For 2 component enthalpy-salinity fields.
 
CoarseAverage m_coarseAverageVector
 For averaging vectors to coarse grids.
 
bool m_timestepReduced
 For reducing the timestep after a failed timestep.
 
int m_dtReduction
 Factor reduce $ \Delta t $ by after a failed timestep.
 
bool m_timestepFailed
 Has this timestep failed?
 
Real m_adv_vel_centering
 If this level has only just been added due to regridding. More...
 
Real m_adv_vel_centering_growth
 Rate at which we grow m_adv_vel_centering each timestep.
 
QuadCFInterp m_quadCFInterpScalar
 For doing quadratic interpolation at coarse fine boundaries (for scalars)
 
QuadCFInterp m_quadCFInterpVector
 For doing quadratic interpolation at coarse fine boundaries (for vectors)
 
PiecewiseLinearFillPatch m_piecewiseLinearFillPatchScalarOne
 For filling one ghost cell around the edge of a set of boxes on a level.
 
PiecewiseLinearFillPatch m_piecewiseLinearFillPatchScalarTwo
 For filling two ghost cells around the edge of a set of boxes on a level.
 
PiecewiseLinearFillPatch m_piecewiseLinearFillPatchScalarThree
 For filling three ghost cells around the edge of a set of boxes on a level.
 
PiecewiseLinearFillPatch m_piecewiseLinearFillPatchScalarFour
 For filling four ghost cells around the edge of a set of boxes on a level.
 
int m_numGhost
 Number of ghost cells for standard calculations.
 
int m_numGhostAdvection
 Number of ghost cells for advection calculations. More...
 
PatchGodunov m_patchGodVelocity
 Handles advection of velocity.
 
PatchGodunov m_patchGodHC
 Operator For multi component H-C advection solve.
 
PatchGodunov m_patchGodTSl
 Operator For multi component T-Sl advection solve.
 
PatchGodunov m_patchGodH
 Operator For single component H advection solve.
 
PatchGodunov m_patchGodC
 Operator For single component C advection solve.
 
PatchGodunov m_patchGodT
 Operator For single component T advection solve.
 
PatchGodunov m_patchGodSl
 Operator For single component Sl advection solve.
 
AdvectionPhysics m_advectionPhysicsVelocity
 Physics for velocity advection.
 
AdvectionPhysics m_advPhysHC
 Physics for coupled H-C advection solve.
 
AdvectionPhysics m_advPhysTSl
 Physics for coupled T-Sl advection solve.
 
AdvectionPhysics m_advPhysH
 Physics for H advection solve.
 
AdvectionPhysics m_advPhysC
 Physics for C advection solve.
 
AdvectionPhysics m_advPhysT
 Physics for T advection solve.
 
AdvectionPhysics m_advPhysSl
 Physics for Sl advection solve.
 
RefCountedPtr< PatchGodunovm_patchGodScalars [m_numScalarVars]
 Handles advection of scalar vars.
 
PhysIBCm_scalarIBC [m_numScalarVars]
 Initial and Boundary Conditions for scalar vars.
 
Projector m_projection
 Handles pressure projection.
 
Projector m_projectionBackup
 Contains the old pressure, in case we need to restart.
 
RefCountedPtr< LevelFluxRegisterm_fluxRegisters [m_numScalarVars]
 Flux registers for scalar variables.
 
RefCountedPtr< LevelFluxRegisterm_fluxRegHC
 Flux Register for multi component HC solves.
 
RefCountedPtr< LevelFluxRegisterm_vectorFluxRegisters [m_numVectorVars]
 Flux registers for vector variables.
 
LevelDomainFluxRegister m_saltDomainFluxRegister
 Flux register for salt along the domain edge.
 
LevelDomainFluxRegister m_heatDomainFluxRegister
 Flux register for heat (enthalpy) along the domain edge.
 
bool m_makeFluxRegForScalarVar [m_numScalarVars]
 Which scalar vars should we bother making flux registers for?
 
bool m_makeFluxRegForVectorVar [m_numVectorVars]
 Which vector vars should we bother making flux registers for?
 
bool m_hasCoarser
 Does a coarser AMRLevelMushyLayer pointer exist?
 
bool m_hasFiner
 Does a finer AMRLevelMushyLayer pointer exist?
 
DisjointBoxLayout m_grids
 Grids on this level.
 
bool m_newGrids_different
 Are the new grids different from current ones?
 
bool m_regrid_smoothing_done
 Has post regrid smoothing been done?
 
Real s_regrid_smoothing_coeff
 Coefficient for post regrid smoothing.
 
bool s_implicit_reflux
 Do reflux implicitly?
 
- Protected Attributes inherited from AMRLevel
ProblemDomain m_problem_domain
 
Vector< Boxm_level_grids
 
int m_level
 
int m_ref_ratio
 
Real m_initial_dt_multiplier
 
Real m_dt
 
Real m_time
 
AMRLevelm_coarser_level_ptr
 
AMRLevelm_finer_level_ptr
 
bool m_isDefined
 

Static Protected Attributes

static BiCGStabSolver< LevelData< FArrayBox > > s_botSolverUStar
 bottom solver for velocity solve multigrid
 
static RelaxSolver< LevelData< FArrayBox > > s_botSolverHC
 Bottom solve for enthalpy-bulk concentration multigrid.
 
- Static Protected Attributes inherited from AMRLevel
static int s_verbosity
 

Additional Inherited Members

- Static Public Member Functions inherited from AMRLevel
static int verbosity ()
 
static void verbosity (int a_verbosity)
 

Detailed Description

AMRLevel for mushy layer calculations.

Big class to manage integration of nonlinear mushy layer equations on an AMR level

Member Enumeration Documentation

◆ solverRestart

What to do when restarting after a failed timestep.

Enumerator
m_restartHalveDt 

Restart and set dt = dt/2.

m_restartResetData 

Restart and return solution to it's state before this timestep.

Member Function Documentation

◆ addHeatSource()

void AMRLevelMushyLayer::addHeatSource ( LevelData< FArrayBox > &  src)
protected

Add a fixed heat source.

$ Q = \frac{Q_0}{\sigma \sqrt{2 \pi}} \exp\left[ - 0.5 \left( \frac{x-x_c}{\sigma} \right)^2 \right] 0.5 \left( 1 + \tanh\left[10 (z-(H-h)) \right]) \right) $

where the parameters are found from the inputs file as: Q_0 is given by heatSource.size sigma is heatSource.width x_c is heatSource.xpos h is heatSource.depth, and measures the depth that the heat source propagates down from the top of the domain

H is the domain height

◆ addMeltPond()

void AMRLevelMushyLayer::addMeltPond ( int  depth,
Real  salinity,
Real  enthalpy,
bool  rescaleExistingSolution = false 
)

Add a melt pond at the top of the domain.

We just do a very simple linear rescaling for now, so the old domain 0 < z < h becomes 0 < z < h-dh, and new variables are denoted by a , i.e. we need a mapping for H(0 < z < h) -> H(0 < z < h-dh), and S->S` We assume that the boundary conditions on the new domain are the same as the old domain at the new extent, i.e. H(h-dh) = H(h), H(h-dh) = S(h). then, H`(z) = H(z*(h-dh)/h) = H(z*(1-dh/h))

Stretching maps z->z, via z = z*stretching

◆ advectActiveTracer()

void AMRLevelMushyLayer::advectActiveTracer ( )
protected

Advect active tracer $ \zeta $.

$ \frac{\partial \zeta}{\partial t} + \mathbf{U} \cdot \nabla (\zeta/\chi) = S(x,z)$ where $ S(x,z) $ is a source/sink of the tracer defined in computeActiveTracerSourceSink()

◆ advectPassiveTracer()

void AMRLevelMushyLayer::advectPassiveTracer ( )
protected

Advect a passive tracer $ \xi $.

$ \frac{\partial \xi}{\partial t} + \mathbf{U} \cdot \nabla (\xi/\chi) = S(x,z)$

◆ advectScalar()

void AMRLevelMushyLayer::advectScalar ( const int  a_scalarVar,
const int  a_advectionVar,
LevelData< FluxBox > &  a_advVel,
bool  doFRupdates = true 
)
protected

Update a scalar field by advection.

a_scalarVar is the field which is updated, whilst a_advectionVar is the field that has been advected. These don't need to be the same, i.e. could do advectScalar(enthalpy, temperature, ...) in order to compute an update for the equation d(enthalpy)/dt + u.grad(temperature) = 0

◆ AMRRefluxLambda()

void AMRLevelMushyLayer::AMRRefluxLambda ( )
protected

Reflux $ \lambda $ across AMR hierarchy.

Only call from base level

◆ calculateTimeIndAdvectionVel()

void AMRLevelMushyLayer::calculateTimeIndAdvectionVel ( Real  time,
LevelData< FluxBox > &  a_advVel 
)
protected

Calculate advection velocity for momentum equations which are time independent.

This source file contains methods relating to the velocity field

◆ computeCCvelocity()

void AMRLevelMushyLayer::computeCCvelocity ( const LevelData< FArrayBox > &  advectionSourceTerm,
Real  a_oldTime,
Real  a_dt,
bool  doFRupdates = false,
bool  a_doProjection = true,
bool  compute_uDelU = true,
bool  a_MACprojection = false 
)
protected

Compute new cell centred velocity.

Generally done via an implicit update due to the viscous and darcy terms. Not required if just considering the Darcy equation.

◆ computeMushDepth()

Real AMRLevelMushyLayer::computeMushDepth ( Real  a_porosity_criteria = 0.99)
protected

Compute mushy layer depth.

Compute horizontally averaged porosity, then determine that the mushy layer starts at the minimum height where the porosity is < a_porosity_criteria.

◆ computeScalarAdvectiveFlux() [1/2]

void AMRLevelMushyLayer::computeScalarAdvectiveFlux ( LevelData< FluxBox > &  a_edgeScal,
LevelData< FArrayBox > &  a_old_scal,
LevelData< FluxBox > &  a_adv_vel,
LevelData< FluxBox > &  a_inflowOutflowVel,
LevelData< FArrayBox > &  a_old_vel,
LevelData< FArrayBox > &  a_diffusiveSrc,
PatchGodunov a_patchGod,
Real  a_old_time,
Real  a_dt 
)
protected

Compute the advective flux given the scalar field, advection velocity and source term specified.

Should work for multiple components

◆ computeScalarAdvectiveFlux() [2/2]

void AMRLevelMushyLayer::computeScalarAdvectiveFlux ( LevelData< FluxBox > &  a_edgeScal,
int  a_advectionVar,
int  a_diffusionVar,
LevelData< FluxBox > &  a_advVel,
Real  a_old_time,
Real  a_dt 
)
protected

Compute the advective flux.

advection var is what we advect diffusion var determines the variable used to calculate the diffusive src. Set to -1 to use no source.

◆ computeScalarAdvectiveFluxMultiComp()

void AMRLevelMushyLayer::computeScalarAdvectiveFluxMultiComp ( LevelData< FluxBox > &  a_edgeScal,
LevelData< FluxBox > &  a_advVel,
PatchGodunov a_patchGod,
LevelData< FArrayBox > &  a_scalOld,
Real  a_old_time,
Real  a_dt 
)
protected

Compute the advective flux for the enthalpy and salinity equations,.

Compute temperature and salinity advective flux.

edgeScal contains two components: $ (\mathbf{U} T, \mathbf{U} S_l) $ or $ (V H, V S) $ depending on the parameters given

◆ computeScalarConcInLiquid()

void AMRLevelMushyLayer::computeScalarConcInLiquid ( LevelData< FArrayBox > &  liquid_tracer_conc,
int  a_tracerVar 
)
protected

Compute concentration of scalar (with bulk concentration $ \xi $ ) in the liquid phase.

$ \xi/\chi $

◆ computeScalarDiffusiveSrc()

void AMRLevelMushyLayer::computeScalarDiffusiveSrc ( int  a_scalarBulkConc,
LevelData< FArrayBox > &  a_src 
)
protected

Compute diffusive src for some scalar with bulk concentration $ \xi $.

$ S = \nabla \cdot \chi \nabla (\xi/\chi)$

◆ computeTotalAdvectiveFluxes()

void AMRLevelMushyLayer::computeTotalAdvectiveFluxes ( LevelData< FluxBox > &  edgeScalTotal)
protected

Compute total fluxes that provide the source term for implicit enthalpy-salinity updates.

$ (\mathbf{U} T + V H \mathbf{k}, \mathbf{U} S_l + V S \mathbf{k})$

◆ computeUstar()

void AMRLevelMushyLayer::computeUstar ( LevelData< FArrayBox > &  a_uDelu,
const LevelData< FArrayBox > &  advectionSourceTerm,
Real  a_oldTime,
Real  a_dt,
bool  doFRupdates,
bool  a_MACprojection = false,
bool  compute_uDelU = true 
)
protected

Compute $ \mathbf{U}^* $, the velocity before we project it.

In particular, we explicitly solve for the velocity $ \left[ 1 - \Delta t^\ell \, \left(Da \, \left(\nabla^\ell\right)^2 - \frac{\chi^{\ell, n+1}}{\Pi^{\ell, n+1}} \right) \right] \mathbf{U}^{\ell, n+1}_* = \mathbf{U}^{\ell, n} + \Delta t^\ell \left\{- \chi^{\ell, n+1/2} \, \nabla^\ell \pi^{\ell, n-1/2} + \chi^{\ell, n+1/2} \left[Pr \, Da^2 \left( Ra_T \theta^{\ell, n+1/2} - Ra_C \Theta^{\ell, n+1/2} \right) \right] - \mathbf{U}_{AD, CC}^{\ell, n+1/2} \cdot \nabla^\ell \left( \frac{\mathbf{U}}{\chi} \right)^{\ell, n+1/2} \right\}$

◆ define() [1/2]

void AMRLevelMushyLayer::define ( MushyLayerOptions  a_opt,
MushyLayerParams  a_params 
)

Defines this AMRLevelMushyLayer.

Porosity for Darcy-Brinkman, permeability for Darcy

◆ define() [2/2]

void AMRLevelMushyLayer::define ( AMRLevel a_coarserLevelPtr,
const ProblemDomain a_problemDomain,
int  a_level,
int  a_refRatio 
)
virtual

Define new AMRLevelMushyLayer from coarser.

Bio

Reimplemented from AMRLevel.

◆ doRegularisationOps()

void AMRLevelMushyLayer::doRegularisationOps ( LevelData< FArrayBox > &  a_scal,
int  a_var,
int  a_comp = 0 
)
protected

Operations to keep variables bounded correctly.

Porosity $ \chi > 0 $ Liquid concentration $ -1 <= \Theta_l <= 0 $

◆ doVelocityAdvection()

bool AMRLevelMushyLayer::doVelocityAdvection ( )
protected

Do we need to do velocity advection?

I.e. do we solve for $ \mathbf{U} \cdot \nabla (\mathbf{U}/\chi) $ terms?

◆ fillBuoyancy()

void AMRLevelMushyLayer::fillBuoyancy ( FArrayBox buoyancy,
FArrayBox temperature,
FArrayBox liquidConc,
FArrayBox porosity,
FArrayBox bodyForce 
)
protected

Compute the buoyancy force.

$ \mathbf{B} = \chi ( Ra_T \theta - Ra_C \Theta_l + \mathbf{F} \cdot \mathbf{k} ) \mathbf{k} $

◆ fillUnprojectedDarcyVelocity()

void AMRLevelMushyLayer::fillUnprojectedDarcyVelocity ( LevelData< FluxBox > &  a_advVel,
Real  time 
)
protected

Fill unprojected darcy velocity.

$ \mathbf{U}^* = (\Pi/\nu) * (Ra_T \theta - Ra_C \Theta_l + \mathbf{F}_B) $

where $ \mathbf{F}_B $ is a body force (=0 by default)

◆ getCoarseScalarDataPointers()

void AMRLevelMushyLayer::getCoarseScalarDataPointers ( const int  a_scalarVar,
LevelData< FArrayBox > **  a_coarserDataOldPtr,
LevelData< FArrayBox > **  a_coarserDataNewPtr,
LevelFluxRegister **  a_coarserFRPtr,
LevelFluxRegister **  a_finerFRPtr,
Real a_tCoarserOld,
Real a_tCoarserNew 
)
protected

Smooth the velocity field.

Utility function to get pointers to coarse level objects

◆ grids()

DisjointBoxLayout AMRLevelMushyLayer::grids ( )

Whether or not to compute diagnostics.

Get the grids on this level

◆ initTimeIndependentPressure()

void AMRLevelMushyLayer::initTimeIndependentPressure ( AMRLevelMushyLayer lev,
int  a_max_num_iter = -1 
)
protected

Initialise phi (pressure for Darcy equation)

By default a_max_num_iter = -1, meaning we use the globally defined value in m_opt.num_init_passes

◆ resetLambda()

void AMRLevelMushyLayer::resetLambda ( )
protected

Flag if a new level has been added.

Reset $ \lambda = 1$

◆ reshapeData()

void AMRLevelMushyLayer::reshapeData ( DisjointBoxLayout  newGrids,
ProblemDomain  newDomain 
)

Transfer data to new disjoint box layout.

Changed to allow changes in either x or y directions.

◆ setDefaults()

void AMRLevelMushyLayer::setDefaults ( )

Set default values.

This source file contains methods for initialising and defining objects

◆ setEta()

void AMRLevelMushyLayer::setEta ( Real  a_eta)

Set eta on all projection operators.

Only call from level 0

◆ smoothScalarField()

void AMRLevelMushyLayer::smoothScalarField ( LevelData< FArrayBox > &  a_phi,
int  a_var,
Real  a_smoothing 
)

Smooth some scalar field.

smoothing done implicitly via an elliptic solve. Obtains smooth field $\phi_{smooth}$ by solving $ (1-\alpha \nabla^2)\phi_{smooth} = \phi $ where $\alpha$ is the smoothing parameter a_smoothing

Member Data Documentation

◆ m_adv_vel_centering

Real AMRLevelMushyLayer::m_adv_vel_centering
protected

If this level has only just been added due to regridding.

The advection velocity is calculated at old_time + m_adv_vel_centering*dt Usually 0.5, but set to something smaller after regridding for stability

◆ m_computedCFL

Real AMRLevelMushyLayer::m_computedCFL
protected

The CFL number we're currently running at.

computed from the current max velocity and the grid spacing.

◆ m_enforceGradP

bool AMRLevelMushyLayer::m_enforceGradP
protected

Should we enforce some analytic $ \nabla P $.

Usually false, but can be useful to do this for benchmark problems

◆ m_maxLambda

Real AMRLevelMushyLayer::m_maxLambda
protected

Store max value of $ \lambda-1 $ from previous timestep.

We may decide to increase $ \eta $ if we're not decreasing lambda quickly enough

◆ m_numGhostAdvection

int AMRLevelMushyLayer::m_numGhostAdvection
protected

Number of ghost cells for advection calculations.

Generally need extra ghost cells when doing advection problems

◆ m_opt

MushyLayerOptions AMRLevelMushyLayer::m_opt
protected

Struct to contain all the options for this simulation.

This is so we can easily pass all the options around, and don't have to keep going into parmparse

◆ m_pressureScaleVar

int AMRLevelMushyLayer::m_pressureScaleVar
protected

Which field scales the pressure.

We either have $ \chi \nabla P $ or $ \Pi \nabla P $

◆ m_usePrevPressureForUStar

bool AMRLevelMushyLayer::m_usePrevPressureForUStar
protected

Should we use previous $ \nabla P $ to calculate $ \mathbf{u}^* $.

We then subtract off $ \nabla P $ before actually doing the projection


The documentation for this class was generated from the following files: