Main Page   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Compound Members   File Members  

AMRNavierStokes Class Reference

#include <AMRNavierStokes.H>

Inheritance diagram for AMRNavierStokes:

Inheritance graph
[legend]
Collaboration diagram for AMRNavierStokes:

Collaboration graph
[legend]
List of all members.

Public Methods

 AMRNavierStokes ()
 AMRNavierStokes (AMRLevel *a_coarser_level_ptr, const Box &a_prob_domain, int a_level, int a_ref_ratio)
 AMRNavierStokes (AMRLevel *a_coarser_level_ptr, const ProblemDomain &a_prob_domain, int a_level, int a_ref_ratio)
virtual ~AMRNavierStokes ()
virtual AMRLevelmakeAMRLevel (AMRLevel *a_coarser_level_ptr, const Box &a_problem_domain, int a_level, int a_ref_ratio) const
virtual AMRLevelmakeAMRLevel (AMRLevel *a_coarser_level_ptr, const ProblemDomain &a_problem_domain, int a_level, int a_ref_ratio) const
virtual void define (AMRLevel *a_coarse_level_ptr, const Box &a_problem_domain, int a_level, int a_ref_ratio)
virtual void define (AMRLevel *a_coarse_level_ptr, const ProblemDomain &a_problem_domain, int a_level, int a_ref_ratio)
virtual Real advance ()
 advance by one timestep (returns what?)

virtual void postTimeStep ()
 things to do after a timestep

virtual void tagCells (IntVectSet &a_tags)
 create tags

virtual void tagCellsInit (IntVectSet &a_tags)
 create tags at initialization

virtual void regrid (const Vector< Box > &a_new_grids)
 regrid

virtual void postRegrid (int a_lBase)
 perform any post-regridding ops -- lBase is the finest unchanged level

virtual void initialGrid (const Vector< Box > &a_new_grids)
 initialize grids

virtual void initialData ()
 initialize data

virtual void postInitialize ()
 things do do after initialization

virtual Real computeDt ()
 compute dt

virtual Real computeInitialDt ()
 compute dt with initial data

virtual void CFL (Real a_cfl)
 set CFL number

virtual void refinementThreshold (Real a_refine_threshold)
 set refinement threshold

bool finestLevel () const
 is this this finest level?

bool isEmpty () const
 is this an empty level?

void computeVorticity (LevelData< FArrayBox > &a_vorticity) const
Real Dx () const
 access functions

Real Nu () const
 viscous coefficient

AMRNavierStokes * finerNSPtr () const
 encapsulates the finer-level-pointer casting

AMRNavierStokes * crseNSPtr () const
 encapsulates the coarser-level pointer casting

LevelData< FArrayBox > & newVel ()
const LevelData< FArrayBox > & newVel () const
LevelData< FArrayBox > & oldVel ()
const LevelData< FArrayBox > & oldVel () const
void velocity (LevelData< FArrayBox > &a_vel, Real a_time) const
 returns velocity at time t

LevelData< FArrayBox > & newLambda ()
const LevelData< FArrayBox > & newLambda () const
LevelData< FArrayBox > & oldLambda ()
const LevelData< FArrayBox > & oldLambda () const
void lambda (LevelData< FArrayBox > &a_lambda, Real a_time) const
 returns lambda at time t

LevelData< FArrayBox > & newScal (const int a_comp)
const LevelData< FArrayBox > & newScal (const int a_comp) const
LevelData< FArrayBox > & oldScal (const int a_comp)
const LevelData< FArrayBox > & oldScal (const int a_comp) const
void scalars (LevelData< FArrayBox > &a_scalars, const Real a_time, const int a_comp) const
 returns scalars at time t

void setPhysBC (const PhysBCUtil &a_BC)
 sets physical BC object

PhysBCUtilgetPhysBCPtr () const
 gets physical BC object


Protected Methods

DisjointBoxLayout loadBalance (const Vector< Box > &a_grids)
void levelSetup (const DisjointBoxLayout &level_domain)
void readParameters ()
 read parameters from parmParse database

void finestLevel (bool a_finest_level)
 set whether this is the finest level or not

void swapOldAndNewStates ()
 move new-time state into old-time state -- also advances time by dt

void resetStates (const Real a_time)
 exchange old and new states, resets time to a_time

void computeAdvectionVelocities (LevelData< FluxBox > &a_adv_vel)
 compute advection velocities

void predictVelocities (LevelData< FArrayBox > &a_uDelU, LevelData< FluxBox > &a_advVel)
void computeUStar (LevelData< FArrayBox > &a_uStar)
void computeUStarInit (LevelData< FArrayBox > &a_uStar)
void advectScalar (LevelData< FArrayBox > &a_new_scal, LevelData< FArrayBox > &a_old_scal, LevelData< FluxBox > &a_adv_vel, LevelFluxRegister *a_crseFluxRegPtr, LevelFluxRegister &a_fineFluxReg, Real a_dt)
 do scalar advection -- assumes all BC's already done

void advectDiffuseScalar (LevelData< FArrayBox > &a_new_scal, LevelData< FArrayBox > &a_old_scal, LevelData< FluxBox > &a_adv_vel, const Real &a_diffusive_coeff, const LevelData< FArrayBox > *a_crse_scal_old, const LevelData< FArrayBox > *a_crse_scal_new, const Real a_old_crse_time, const Real a_new_crse_time, LevelFluxRegister *a_crseFluxRegPtr, LevelFluxRegister &a_fineFluxReg, DomainGhostBC &a_scalPhysBC, Real a_dt)
 advect/diffuse scalar -- assumes all BC's already done?

void initializeGlobalPressure ()
 manages pressure initialization after init or regrid

void initializeLevelPressure (Real a_currentTime, Real a_dtInit)
 does level-based pressure initalization

void smoothVelocityField (int lbase)
 smooths composite velocity field after regridding

void computeLapVel (LevelData< FArrayBox > &a_lapVel, LevelData< FArrayBox > &a_vel, const LevelData< FArrayBox > *a_crseVelPtr)
 compute Laplacian(velocity)

void computeLapScal (LevelData< FArrayBox > &a_lapScal, LevelData< FArrayBox > &a_scal, DomainGhostBC &a_physBC, const LevelData< FArrayBox > *a_crseScalPtr)
 compute Laplacian(scalars)

void fillVelocity (LevelData< FArrayBox > &a_vel, Real a_time)
 fill grown velocity field using piecewise-constant interp

void fillVelocity (LevelData< FArrayBox > &a_vel, Real a_time, int vel_comp, int dest_comp, int num_comp)
 fill grown velocity field using piecewise-constant interp

void fillLambda (LevelData< FArrayBox > &a_lambda, Real a_time) const
void fillScalars (LevelData< FArrayBox > &a_scal, Real a_time, const int a_comp) const
 fill interior and boundary cells (assume a_scal already sized)

virtual int numPlotComps () const
 number of components in plotfiles

virtual void getPlotData (LevelData< FArrayBox > &a_plot_data) const
 stuff plotfile data

void defineRegridAMRSolver (AMRSolver &a_solver, const Vector< DisjointBoxLayout > &a_grids, const Vector< ProblemDomain > &a_domains, const Vector< Real > &a_amrDx, const Vector< int > &a_refRatios, const int &a_lBase)
 defines (helmholtz) AMRSolver used in smoothing during regridding


Protected Attributes

LevelData< FArrayBox > * m_vel_old_ptr
 velocity vector at old time

LevelData< FArrayBox > * m_vel_new_ptr
 velocity vector at new time

LevelData< FArrayBox > * m_lambda_old_ptr
 freestream preservation qty at old time

LevelData< FArrayBox > * m_lambda_new_ptr
 freestream preservation qty at new time

Vector< LevelData< FArrayBox > * > m_scal_old
 scalars at old time

Vector< LevelData< FArrayBox > * > m_scal_new
 scalars at new time

Vector< Tuple< Real, SpaceDim > > m_dataPoints
 points where data is dumped out

Real m_dx
 grid spacing

Real m_cfl
 cfl number

Real m_dt_save
 dt from most recent timestep

int m_level_steps
 number of steps on this level

CoarseAverage m_coarse_average
 averaging from fine to coarse level for velocity

CoarseAverage m_coarse_average_lambda
 averaging from fine to coarse level for freestream preservation qty

CCProjector m_projection
 cell-centered projection object

LevelFluxRegister m_flux_register
 momentum flux register

LevelFluxRegister m_lambda_flux_reg
 freestream preservation qty flux register

Vector< LevelFluxRegister * > m_scal_fluxreg_ptrs
 conserved scalar flux register

QuadCFInterp m_velCFInterp
 does quadratic coarse-fine interpolation for velocities

VelocityPoissonOp m_velocityPoissonOp
 PoissonOp for velocities.

PoissonOp m_scalarsPoissonOp
 PoissonOp for advected scalars (for diffusion).

LevelHelmholtzSolver m_viscous_solver
 levelsolver for viscous solves

Real m_refine_threshold
 refinement threshold for regridding

bool m_finest_level
 true if this is the finest extant level

bool m_is_empty
 true if this level is an empty level (no grids)

bool m_regrid_smoothing_done
 internal indicator whether post-regrid smoothing stuff has been done

PhysBCUtilm_physBCPtr
 physical BC object


Static Protected Attributes

const int s_num_vel_comps = SpaceDim
 number of components of m_state (vel comps + lambda)

int s_num_scal_comps
 number of scalars

const int s_max_scal_comps = SpaceDim
 maximum possible number of scalar components

const char * s_vel_names [s_num_vel_comps]
 names of components

const char * s_scal_names [s_max_scal_comps]
 names of scalars

Vector< Reals_scal_coeffs
 diffusion coefficients for scalar advection-diffusion

Vector< Reals_domLength
 size of physical domain

bool s_ppInit
 have we initialized static variables from parmParse?

Real s_init_shrink
 factor to shrink initial timestep

Real s_max_dt
 maximum allowable timestep

Real s_prescribedDt
bool s_project_initial_vel
 should we project initial velocity field (default is yes)

bool s_initialize_pressures
 should we initialize pressures after grid generation & regridding?

bool s_smooth_after_regrid
 do antidiffusive smoothing after regridding?

bool s_dump_smoothing_plots
 dump before-and-after plotfiles for post-regrid smoothing (for debug)

Real s_regrid_smoothing_coeff
 antidiffusive smoothing coefficient (multiplies nu*dtCrse)

int s_num_init_passes
 number of times to iterate when initializing pressures

bool s_reflux_momentum
 should we reflux momentum?

bool s_set_bogus_values
 debugging option -- initially set arrays to s_bogus_value

Real s_bogus_value
 debugging option -- value to which arrays are initially set

bool s_tag_vorticity
 tag on vorticity?

bool s_tag_theta
 tag on scalar (theta)?

Real s_vort_factor
 factor for tagging -- tag if vorticity greater than factor*max_vort

int s_tags_grow
 amount by which to grow tagged regions before calling meshrefine

bool s_implicit_reflux
 do implicit refluxing for velocity?

bool s_applyFreestreamCorrection
 apply freestream preservation correction?

bool s_reflux_scal
 should we reflux scalars?

bool s_implicit_scal_reflux
 do implicit refluxing for advected/diffused scalars?

int s_viscous_solver_type
 viscous solver type: 0 = backwards euler, 1=Crank-Nicolson, 2=TGA

Real s_viscous_solver_tol
 viscous solver tolerance (defaults to 1e-7)

Real s_nu
 viscosity

bool s_boussinesq
 use Boussinesq approximation -- theta (temperature) provides forcing

Real s_theta_coeff
 factor to use in theta forcing

Real s_rayleighNo
 Rayleigh number, if boussinesq.

Real s_prandtlNo
 Prandtl number, if boussinesq.

Real s_wallTemp
 wall temperature, if boussinesq

bool s_specifyInitialGrids
 do we specify initial grids in a gridFile?

string s_initialGridFile
 gridfile where initial grids are spec'd, if relevant

bool s_init_vel_from_vorticity
 initialize velocity field from vorticity field

Real s_backgroundVel
 background velocity (in the z-direction)


Constructor & Destructor Documentation

AMRNavierStokes::AMRNavierStokes  
 

AMRNavierStokes::AMRNavierStokes AMRLevel   a_coarser_level_ptr,
const Box   a_prob_domain,
int    a_level,
int    a_ref_ratio
 

AMRNavierStokes::AMRNavierStokes AMRLevel   a_coarser_level_ptr,
const ProblemDomain   a_prob_domain,
int    a_level,
int    a_ref_ratio
 

virtual AMRNavierStokes::~AMRNavierStokes   [virtual]
 


Member Function Documentation

virtual Real AMRNavierStokes::advance   [virtual]
 

advance by one timestep (returns what?)

void AMRNavierStokes::advectDiffuseScalar LevelData< FArrayBox > &    a_new_scal,
LevelData< FArrayBox > &    a_old_scal,
LevelData< FluxBox > &    a_adv_vel,
const Real   a_diffusive_coeff,
const LevelData< FArrayBox > *    a_crse_scal_old,
const LevelData< FArrayBox > *    a_crse_scal_new,
const Real    a_old_crse_time,
const Real    a_new_crse_time,
LevelFluxRegister   a_crseFluxRegPtr,
LevelFluxRegister   a_fineFluxReg,
DomainGhostBC   a_scalPhysBC,
Real    a_dt
[protected]
 

advect/diffuse scalar -- assumes all BC's already done?

void AMRNavierStokes::advectScalar LevelData< FArrayBox > &    a_new_scal,
LevelData< FArrayBox > &    a_old_scal,
LevelData< FluxBox > &    a_adv_vel,
LevelFluxRegister   a_crseFluxRegPtr,
LevelFluxRegister   a_fineFluxReg,
Real    a_dt
[protected]
 

do scalar advection -- assumes all BC's already done

virtual void AMRNavierStokes::CFL Real    a_cfl [virtual]
 

set CFL number

void AMRNavierStokes::computeAdvectionVelocities LevelData< FluxBox > &    a_adv_vel [protected]
 

compute advection velocities

virtual Real AMRNavierStokes::computeDt   [virtual]
 

compute dt

Implements AMRLevel.

virtual Real AMRNavierStokes::computeInitialDt   [virtual]
 

compute dt with initial data

Implements AMRLevel.

void AMRNavierStokes::computeLapScal LevelData< FArrayBox > &    a_lapScal,
LevelData< FArrayBox > &    a_scal,
DomainGhostBC   a_physBC,
const LevelData< FArrayBox > *    a_crseScalPtr
[protected]
 

compute Laplacian(scalars)

void AMRNavierStokes::computeLapVel LevelData< FArrayBox > &    a_lapVel,
LevelData< FArrayBox > &    a_vel,
const LevelData< FArrayBox > *    a_crseVelPtr
[protected]
 

compute Laplacian(velocity)

void AMRNavierStokes::computeUStar LevelData< FArrayBox > &    a_uStar [protected]
 

void AMRNavierStokes::computeUStarInit LevelData< FArrayBox > &    a_uStar [protected]
 

void AMRNavierStokes::computeVorticity LevelData< FArrayBox > &    a_vorticity const
 

AMRNavierStokes* AMRNavierStokes::crseNSPtr   const
 

encapsulates the coarser-level pointer casting

virtual void AMRNavierStokes::define AMRLevel   a_coarse_level_ptr,
const ProblemDomain   a_problem_domain,
int    a_level,
int    a_ref_ratio
[virtual]
 

Defines this AMRLevel.

{\bf Arguments:}\ a_coarser_level_ptr (not modified): pointer to next coarser level object.\ a_problem_domain (not modified): problem domain of this level.\ a_level (not modified): index of this level. The base level is zero.\ a_ref_ratio (not modified): the refinement ratio between this level and the next finer level.

Reimplemented from AMRLevel.

virtual void AMRNavierStokes::define AMRLevel   a_coarse_level_ptr,
const Box   a_problem_domain,
int    a_level,
int    a_ref_ratio
[virtual]
 

Defines this AMRLevel.

{\bf Arguments:}\ a_coarser_level_ptr (not modified): pointer to next coarser level object.\ a_problem_domain (not modified): problem domain of this level.\ a_level (not modified): index of this level. The base level is zero.\ a_ref_ratio (not modified): the refinement ratio between this level and the next finer level.

Reimplemented from AMRLevel.

void AMRNavierStokes::defineRegridAMRSolver AMRSolver   a_solver,
const Vector< DisjointBoxLayout > &    a_grids,
const Vector< ProblemDomain > &    a_domains,
const Vector< Real > &    a_amrDx,
const Vector< int > &    a_refRatios,
const int &    a_lBase
[protected]
 

defines (helmholtz) AMRSolver used in smoothing during regridding

Real AMRNavierStokes::Dx   const
 

access functions

void AMRNavierStokes::fillLambda LevelData< FArrayBox > &    a_lambda,
Real    a_time
const [protected]
 

resize temporary lambda for advection tracing computation and fill interior and boundary cells

void AMRNavierStokes::fillScalars LevelData< FArrayBox > &    a_scal,
Real    a_time,
const int    a_comp
const [protected]
 

fill interior and boundary cells (assume a_scal already sized)

void AMRNavierStokes::fillVelocity LevelData< FArrayBox > &    a_vel,
Real    a_time,
int    vel_comp,
int    dest_comp,
int    num_comp
[protected]
 

fill grown velocity field using piecewise-constant interp

void AMRNavierStokes::fillVelocity LevelData< FArrayBox > &    a_vel,
Real    a_time
[protected]
 

fill grown velocity field using piecewise-constant interp

AMRNavierStokes* AMRNavierStokes::finerNSPtr   const
 

encapsulates the finer-level-pointer casting

void AMRNavierStokes::finestLevel bool    a_finest_level [protected]
 

set whether this is the finest level or not

bool AMRNavierStokes::finestLevel   const
 

is this this finest level?

PhysBCUtil* AMRNavierStokes::getPhysBCPtr   const
 

gets physical BC object

virtual void AMRNavierStokes::getPlotData LevelData< FArrayBox > &    a_plot_data const [protected, virtual]
 

stuff plotfile data

virtual void AMRNavierStokes::initialData   [virtual]
 

initialize data

Implements AMRLevel.

virtual void AMRNavierStokes::initialGrid const Vector< Box > &    a_new_grids [virtual]
 

initialize grids

Implements AMRLevel.

void AMRNavierStokes::initializeGlobalPressure   [protected]
 

manages pressure initialization after init or regrid

void AMRNavierStokes::initializeLevelPressure Real    a_currentTime,
Real    a_dtInit
[protected]
 

does level-based pressure initalization

bool AMRNavierStokes::isEmpty   const
 

is this an empty level?

void AMRNavierStokes::lambda LevelData< FArrayBox > &    a_lambda,
Real    a_time
const
 

returns lambda at time t

void AMRNavierStokes::levelSetup const DisjointBoxLayout   level_domain [protected]
 

DisjointBoxLayout AMRNavierStokes::loadBalance const Vector< Box > &    a_grids [protected]
 

virtual AMRLevel* AMRNavierStokes::makeAMRLevel AMRLevel   a_coarser_level_ptr,
const ProblemDomain   a_problem_domain,
int    a_level,
int    a_ref_ratio
const [virtual]
 

virtual AMRLevel* AMRNavierStokes::makeAMRLevel AMRLevel   a_coarser_level_ptr,
const Box   a_problem_domain,
int    a_level,
int    a_ref_ratio
const [virtual]
 

const LevelData<FArrayBox>& AMRNavierStokes::newLambda   const
 

LevelData<FArrayBox>& AMRNavierStokes::newLambda  
 

const LevelData<FArrayBox>& AMRNavierStokes::newScal const int    a_comp const
 

LevelData<FArrayBox>& AMRNavierStokes::newScal const int    a_comp
 

const LevelData<FArrayBox>& AMRNavierStokes::newVel   const
 

LevelData<FArrayBox>& AMRNavierStokes::newVel  
 

Real AMRNavierStokes::Nu   const
 

viscous coefficient

virtual int AMRNavierStokes::numPlotComps   const [protected, virtual]
 

number of components in plotfiles

const LevelData<FArrayBox>& AMRNavierStokes::oldLambda   const
 

LevelData<FArrayBox>& AMRNavierStokes::oldLambda  
 

const LevelData<FArrayBox>& AMRNavierStokes::oldScal const int    a_comp const
 

LevelData<FArrayBox>& AMRNavierStokes::oldScal const int    a_comp
 

const LevelData<FArrayBox>& AMRNavierStokes::oldVel   const
 

LevelData<FArrayBox>& AMRNavierStokes::oldVel  
 

virtual void AMRNavierStokes::postInitialize   [virtual]
 

things do do after initialization

Implements AMRLevel.

virtual void AMRNavierStokes::postRegrid int    a_lBase [virtual]
 

perform any post-regridding ops -- lBase is the finest unchanged level

Reimplemented from AMRLevel.

virtual void AMRNavierStokes::postTimeStep   [virtual]
 

things to do after a timestep

Implements AMRLevel.

void AMRNavierStokes::predictVelocities LevelData< FArrayBox > &    a_uDelU,
LevelData< FluxBox > &    a_advVel
[protected]
 

void AMRNavierStokes::readParameters   [protected]
 

read parameters from parmParse database

virtual void AMRNavierStokes::refinementThreshold Real    a_refine_threshold [virtual]
 

set refinement threshold

virtual void AMRNavierStokes::regrid const Vector< Box > &    a_new_grids [virtual]
 

regrid

Implements AMRLevel.

void AMRNavierStokes::resetStates const Real    a_time [protected]
 

exchange old and new states, resets time to a_time

void AMRNavierStokes::scalars LevelData< FArrayBox > &    a_scalars,
const Real    a_time,
const int    a_comp
const
 

returns scalars at time t

void AMRNavierStokes::setPhysBC const PhysBCUtil   a_BC
 

sets physical BC object

void AMRNavierStokes::smoothVelocityField int    lbase [protected]
 

smooths composite velocity field after regridding

void AMRNavierStokes::swapOldAndNewStates   [protected]
 

move new-time state into old-time state -- also advances time by dt

virtual void AMRNavierStokes::tagCells IntVectSet   a_tags [virtual]
 

create tags

Implements AMRLevel.

virtual void AMRNavierStokes::tagCellsInit IntVectSet   a_tags [virtual]
 

create tags at initialization

Implements AMRLevel.

void AMRNavierStokes::velocity LevelData< FArrayBox > &    a_vel,
Real    a_time
const
 

returns velocity at time t


Member Data Documentation

Real AMRNavierStokes::m_cfl [protected]
 

cfl number

CoarseAverage AMRNavierStokes::m_coarse_average [protected]
 

averaging from fine to coarse level for velocity

CoarseAverage AMRNavierStokes::m_coarse_average_lambda [protected]
 

averaging from fine to coarse level for freestream preservation qty

Vector<Tuple<Real, SpaceDim> > AMRNavierStokes::m_dataPoints [protected]
 

points where data is dumped out

Real AMRNavierStokes::m_dt_save [protected]
 

dt from most recent timestep

Real AMRNavierStokes::m_dx [protected]
 

grid spacing

bool AMRNavierStokes::m_finest_level [protected]
 

true if this is the finest extant level

LevelFluxRegister AMRNavierStokes::m_flux_register [protected]
 

momentum flux register

bool AMRNavierStokes::m_is_empty [protected]
 

true if this level is an empty level (no grids)

LevelFluxRegister AMRNavierStokes::m_lambda_flux_reg [protected]
 

freestream preservation qty flux register

LevelData<FArrayBox>* AMRNavierStokes::m_lambda_new_ptr [protected]
 

freestream preservation qty at new time

LevelData<FArrayBox>* AMRNavierStokes::m_lambda_old_ptr [protected]
 

freestream preservation qty at old time

int AMRNavierStokes::m_level_steps [protected]
 

number of steps on this level

PhysBCUtil* AMRNavierStokes::m_physBCPtr [protected]
 

physical BC object

CCProjector AMRNavierStokes::m_projection [protected]
 

cell-centered projection object

Real AMRNavierStokes::m_refine_threshold [protected]
 

refinement threshold for regridding

bool AMRNavierStokes::m_regrid_smoothing_done [protected]
 

internal indicator whether post-regrid smoothing stuff has been done

Vector<LevelFluxRegister*> AMRNavierStokes::m_scal_fluxreg_ptrs [protected]
 

conserved scalar flux register

Vector<LevelData<FArrayBox>*> AMRNavierStokes::m_scal_new [protected]
 

scalars at new time

Vector<LevelData<FArrayBox>*> AMRNavierStokes::m_scal_old [protected]
 

scalars at old time

PoissonOp AMRNavierStokes::m_scalarsPoissonOp [protected]
 

PoissonOp for advected scalars (for diffusion).

LevelData<FArrayBox>* AMRNavierStokes::m_vel_new_ptr [protected]
 

velocity vector at new time

LevelData<FArrayBox>* AMRNavierStokes::m_vel_old_ptr [protected]
 

velocity vector at old time

QuadCFInterp AMRNavierStokes::m_velCFInterp [protected]
 

does quadratic coarse-fine interpolation for velocities

VelocityPoissonOp AMRNavierStokes::m_velocityPoissonOp [protected]
 

PoissonOp for velocities.

LevelHelmholtzSolver AMRNavierStokes::m_viscous_solver [protected]
 

levelsolver for viscous solves

bool AMRNavierStokes::s_applyFreestreamCorrection [static, protected]
 

apply freestream preservation correction?

Real AMRNavierStokes::s_backgroundVel [static, protected]
 

background velocity (in the z-direction)

Real AMRNavierStokes::s_bogus_value [static, protected]
 

debugging option -- value to which arrays are initially set

bool AMRNavierStokes::s_boussinesq [static, protected]
 

use Boussinesq approximation -- theta (temperature) provides forcing

Vector<Real> AMRNavierStokes::s_domLength [static, protected]
 

size of physical domain

bool AMRNavierStokes::s_dump_smoothing_plots [static, protected]
 

dump before-and-after plotfiles for post-regrid smoothing (for debug)

bool AMRNavierStokes::s_implicit_reflux [static, protected]
 

do implicit refluxing for velocity?

bool AMRNavierStokes::s_implicit_scal_reflux [static, protected]
 

do implicit refluxing for advected/diffused scalars?

Real AMRNavierStokes::s_init_shrink [static, protected]
 

factor to shrink initial timestep

bool AMRNavierStokes::s_init_vel_from_vorticity [static, protected]
 

initialize velocity field from vorticity field

string AMRNavierStokes::s_initialGridFile [static, protected]
 

gridfile where initial grids are spec'd, if relevant

bool AMRNavierStokes::s_initialize_pressures [static, protected]
 

should we initialize pressures after grid generation & regridding?

Real AMRNavierStokes::s_max_dt [static, protected]
 

maximum allowable timestep

const int AMRNavierStokes::s_max_scal_comps = SpaceDim [static, protected]
 

maximum possible number of scalar components

Real AMRNavierStokes::s_nu [static, protected]
 

viscosity

int AMRNavierStokes::s_num_init_passes [static, protected]
 

number of times to iterate when initializing pressures

int AMRNavierStokes::s_num_scal_comps [static, protected]
 

number of scalars

const int AMRNavierStokes::s_num_vel_comps = SpaceDim [static, protected]
 

number of components of m_state (vel comps + lambda)

bool AMRNavierStokes::s_ppInit [static, protected]
 

have we initialized static variables from parmParse?

Real AMRNavierStokes::s_prandtlNo [static, protected]
 

Prandtl number, if boussinesq.

Real AMRNavierStokes::s_prescribedDt [static, protected]
 

bool AMRNavierStokes::s_project_initial_vel [static, protected]
 

should we project initial velocity field (default is yes)

Real AMRNavierStokes::s_rayleighNo [static, protected]
 

Rayleigh number, if boussinesq.

bool AMRNavierStokes::s_reflux_momentum [static, protected]
 

should we reflux momentum?

bool AMRNavierStokes::s_reflux_scal [static, protected]
 

should we reflux scalars?

Real AMRNavierStokes::s_regrid_smoothing_coeff [static, protected]
 

antidiffusive smoothing coefficient (multiplies nu*dtCrse)

Vector<Real> AMRNavierStokes::s_scal_coeffs [static, protected]
 

diffusion coefficients for scalar advection-diffusion

const char* AMRNavierStokes::s_scal_names[s_max_scal_comps] [static, protected]
 

names of scalars

bool AMRNavierStokes::s_set_bogus_values [static, protected]
 

debugging option -- initially set arrays to s_bogus_value

bool AMRNavierStokes::s_smooth_after_regrid [static, protected]
 

do antidiffusive smoothing after regridding?

bool AMRNavierStokes::s_specifyInitialGrids [static, protected]
 

do we specify initial grids in a gridFile?

bool AMRNavierStokes::s_tag_theta [static, protected]
 

tag on scalar (theta)?

bool AMRNavierStokes::s_tag_vorticity [static, protected]
 

tag on vorticity?

int AMRNavierStokes::s_tags_grow [static, protected]
 

amount by which to grow tagged regions before calling meshrefine

Real AMRNavierStokes::s_theta_coeff [static, protected]
 

factor to use in theta forcing

const char* AMRNavierStokes::s_vel_names[s_num_vel_comps] [static, protected]
 

names of components

Real AMRNavierStokes::s_viscous_solver_tol [static, protected]
 

viscous solver tolerance (defaults to 1e-7)

int AMRNavierStokes::s_viscous_solver_type [static, protected]
 

viscous solver type: 0 = backwards euler, 1=Crank-Nicolson, 2=TGA

Real AMRNavierStokes::s_vort_factor [static, protected]
 

factor for tagging -- tag if vorticity greater than factor*max_vort

Real AMRNavierStokes::s_wallTemp [static, protected]
 

wall temperature, if boussinesq


The documentation for this class was generated from the following file:
Generated on Thu Aug 29 11:07:35 2002 for Chombo&INS by doxygen1.2.16