IBAMR  IBAMR version 0.19.
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
SAMRAI::appu::EmbeddedBoundaryGeometry< DIM > Class Template Reference

Class EmbeddedBoundaryGeometry provides embedded boundary mesh construction, storage, and management on an AMR hierarchy.
More...

#include <EmbeddedBoundaryGeometry.h>

Inheritance diagram for SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >:
Inheritance graph
[legend]

Public Types

enum  CELL_TYPE { SOLID = EmbeddedBoundaryDefines::SOLID, CUT = EmbeddedBoundaryDefines::CUT, BORDER = EmbeddedBoundaryDefines::BORDER, FLOW = EmbeddedBoundaryDefines::FLOW }
 
enum  NODE_TYPE { OUTSIDE = EmbeddedBoundaryDefines::OUTSIDE, INSIDE = EmbeddedBoundaryDefines::INSIDE, BOUNDARY = EmbeddedBoundaryDefines::BOUNDARY, ONBOUNDARY = EmbeddedBoundaryDefines::ONBOUNDARY }
 
enum  PACK_RETURN_TYPE { VISIT_ALLZERO = 0, VISIT_ALLONE = 1, VISIT_MIXED = 2 }
 Enumerated type for the allowable return values for packMaterialFractionsIntoDoubleBuffer() and packSpeciesFractionsIntoDoubleBuffer(). More...
 

Public Member Functions

 EmbeddedBoundaryGeometry (const std::string &object_name, tbox::Pointer< tbox::Database > input_db=NULL, const tbox::Pointer< geom::CartesianGridGeometry< DIM > > grid_geom=NULL, const hier::IntVector< DIM > &nghosts=hier::IntVector< DIM >(0))
 
 ~EmbeddedBoundaryGeometry ()
 
void buildEmbeddedBoundaryOnLevel (const tbox::Pointer< hier::PatchLevel< DIM > > level, const tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy=NULL, const tbox::Pointer< hier::PatchLevel< DIM > > old_level=NULL)
 
void tagInsideOutsideNodesOnLevel (const tbox::Pointer< hier::PatchLevel< DIM > > level)
 
void registerVisItDataWriter (tbox::Pointer< appu::VisItDataWriter< DIM > > visit_writer)
 
int getCellFlagDataId () const
 
int getCellVolumeDataId () const
 
int getIndexCutCellDataId () const
 
int getNodeInsideOutsideDataId () const
 
double computeTotalVolumeOnLevel (const tbox::Pointer< hier::PatchLevel< DIM > > level)
 
void setGridGeometry (const tbox::Pointer< geom::CartesianGridGeometry< DIM > > grid_geom)
 
void writeLevelEmbeddedBoundaryDataToFile (const tbox::Pointer< hier::PatchLevel< DIM > > level, const std::string &dirname) const
 
int packMaterialFractionsIntoDoubleBuffer (double *dbuffer, const hier::Patch< DIM > &patch, const hier::Box< DIM > &region, const std::string &material_name) const
 
void setPhysicalBoundaryConditions (hier::Patch< DIM > &patch, const double fill_time, const hier::IntVector< DIM > &ghost_width_to_fill)
 
virtual void preprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio)
 
virtual void postprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio)
 
virtual hier::IntVector< DIM > getRefineOpStencilWidth () const
 
void putToDatabase (tbox::Pointer< tbox::Database > db)
 
virtual int packMaterialFractionsIntoSparseBuffers (int *mat_list, std::vector< int > &mix_zones, std::vector< int > &mix_mat, std::vector< double > &vol_fracs, std::vector< int > &next_mat, const hier::Patch< DIM > &patch, const hier::Box< DIM > &region) const
 Pack sparse volume fraction data. More...
 
virtual int packSpeciesFractionsIntoDoubleBuffer (double *buffer, const hier::Patch< DIM > &patch, const hier::Box< DIM > &region, const std::string &material_name, const std::string &species_name) const
 This function packs cell-centered species fractions for the given species. More...
 
virtual void preprocessRefineBoxes (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::BoxList< DIM > &fine_boxes, const hier::IntVector< DIM > &ratio)
 
virtual void postprocessRefineBoxes (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::BoxList< DIM > &fine_boxes, const hier::IntVector< DIM > &ratio)
 

Private Member Functions

void initializeVariables (const hier::IntVector< DIM > &nghosts)
 
void computeEmbeddedBoundaryOnLevelWithPackage (const tbox::Pointer< hier::PatchLevel< DIM > > level, int &cut_cells_on_level)
 
void computeEmbeddedBoundaryOnLevel (const tbox::Pointer< hier::PatchLevel< DIM > > level, int &cut_cells_on_level, double &l2_volume_error, double &l2_area_error, double &max_volume_error, double &max_area_error)
 
void calculateCutCellInformation (appu::CutCell< DIM > &cut_cell, const double *lower, const double *upper, const double &fullcellvol, const double *fullareas, double &volume_error_estimate, double &area_error_estimate)
 
void calculateBoundaryNodeInformation (appu::CutCell< DIM > &cut_cell, tbox::Pointer< pdat::NodeData< DIM, int > > &node_flag, const double *lower, const double *upper, const double *dx)
 
double calculateVolume (const double *cell_lower, const double *cell_upper, double &error_estimate) const
 
int recursiveCalculateVolume (const double *cell_lower, const double *cell_upper, double &volume, double &error_estimate, int &subdivide_level) const
 
double calculateArea (const double *cell_lower, const double *cell_upper, const int face, double &error_estimate) const
 
int recursiveCalculateArea (const double *cell_lower, const double *cell_upper, const int face, double &area, double &error_estimate, int &subdivide_level) const
 
int classifyCell (const double *lower, const double *upper) const
 
void refineEmbeddedBoundary (const tbox::Pointer< hier::PatchLevel< DIM > > level, const tbox::Pointer< hier::PatchLevel< DIM > > old_level=NULL, const tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy=NULL)
 
void setEmbeddedBoundaryAtPhysicalBoundaries (const tbox::Pointer< hier::PatchLevel< DIM > > level)
 
void setLevelRatioInformation (const tbox::Pointer< hier::PatchLevel< DIM > > level)
 
void doNativeShapeInsideOutside (const int *nx, const double *dx, const double *origin, int *node_flag) const
 
void readLevelEmbeddedBoundaryDataFromFile (const tbox::Pointer< hier::PatchLevel< DIM > > level, const std::string &dirname) const
 
void getFromInput (tbox::Pointer< tbox::Database > db, const bool is_from_restart)
 
void getFromRestart ()
 

Private Attributes

std::string d_object_name
 
tbox::Pointer< geom::CartesianGridGeometry< DIM > > d_grid_geometry
 
tbox::Pointer< appu::VisItDataWriter< DIM > > d_visit_writer
 
bool d_verbose
 
std::string d_flow_type
 
hier::IntVector< DIM > d_ebdry_nghosts
 
tbox::Pointer< pdat::IndexVariable< DIM, appu::CutCell< DIM >, pdat::CellGeometry< DIM > > > d_ebdry_var
 
tbox::Pointer< pdat::CellVariable< DIM, int > > d_cell_flag_var
 
tbox::Pointer< pdat::CellVariable< DIM, double > > d_cell_vol_var
 
tbox::Pointer< pdat::NodeVariable< DIM, int > > d_node_flag_var
 
int d_ebdry_data_id
 
int d_ebdry_scratch_id
 
int d_cell_flag_data_id
 
int d_cell_flag_scratch_id
 
int d_cell_vol_data_id
 
int d_cell_vol_scratch_id
 
int d_node_flag_data_id
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_ebdry_refine_alg
 
int d_max_levels
 
tbox::Array< hier::IntVector< DIM > > d_ratio_to_coarser
 
tbox::Array< hier::IntVector< DIM > > d_level_ratios
 
bool d_use_cubes
 
tbox::Pointer< appu::CubesPatchInterface< DIM > > d_cubes_interface
 
bool d_use_eleven_inside_outside
 
bool d_use_eleven_boundary_node
 
tbox::Pointer< ElevenPatchInterface< DIM > > d_eleven_interface
 
tbox::Array< tbox::Pointer< appu::EmbeddedBoundaryShape< DIM > > > d_shapes
 
int d_max_subdivides
 
bool d_read_ebdry
 
bool d_write_ebdry
 
std::string d_ebdry_dirname
 
bool d_use_recursive_algs
 
bool d_compute_areas_and_normal
 
bool d_compute_cutcell_index_data
 
bool d_compute_boundary_node_data
 
tbox::Array< intd_node_bdry_cond
 
tbox::Array< intd_edge_bdry_cond
 
tbox::Array< intd_face_bdry_cond
 
tbox::Pointer< tbox::Timert_compute_eb
 
tbox::Pointer< tbox::Timert_calc_node_inout
 
tbox::Pointer< tbox::Timert_calc_volume
 
tbox::Pointer< tbox::Timert_calc_area
 
tbox::Pointer< tbox::Timert_calc_boundary_node
 
tbox::Pointer< tbox::Timert_eb_phys_bdry
 
tbox::Pointer< tbox::Timert_read_write_eb
 
int d_step_count
 

Detailed Description

template<int DIM>
class SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >

The embedded boundary may be constructed from a set of analytic shapes supplied through input. The following outlines the steps required:

  1. Construct an EmbeddedBoundaryGeometry object, supplying the input file that contains the shape entries, a pointer to the Cartesian grid geometry, and the desired number of ghosts to be used in defining the embedded boundary:
*    EmbeddedBoundaryGeometry* eb_geom = 
*       new EmbeddedBoundaryGeometry("EmbeddedBoundaryGeometry",
*                                     input_db->getDatabase("EBdryGeometry"),
*                                     grid_geometry,
*                                     nghosts);
*    

Note: The cart_grid_geometry argument is optional and may be supplied later using the "setGridGeometry()" method. However, it must be set before the "buildEmbeddedBoundaryOnLevel()" is called. The nghosts argument is also optional and is set to zero by default.

  1. Build the embedded boundary on the levels of the hierarchy:
*    Pointer<PatchLevel<DIM> > level = hierarchy->getPatchLevel(ln);
*    eb_geom->buildEmbeddedBoundaryOnLevel(level);
*    

Note: It is also possible to pass in a hierarchy and an old level as arguments. If supplied, this information can be used to accelerate construction of the boundary on the supplied level.

  1. Access information about the embedded boundary from patches:
    • The "cell flag" identifies whether a cell is cut, solid, or flow
    • The "node flag" identifies nodes as inside, outside, boundary (first node just inside solid boundary), or on-boundary (within a specified distance from the solid boundary).
    • The "volume fraction" will be 0.0 for solid cells, 1.0 for flow cells, and somewhere in-between for cut cells
    • The list of "cut cells" holds information about the normal, face areas, etc. on specific cells that are cut.
*    int cell_flag_index = eb_geom->getCellFlagDataId();  
*    int node_flag_index = eb_geom->getNodeInsideOutsideDataId();  
*    int vol_frac_index = eb_geom->getCellVolumeDataId();  
*    int cut_cell_index = eb_geom->getIndexCutCellDataId();
*
*    tbox::Pointer<CellData<DIM,int> > cell_flag_data = 
*        patch->getPatchData(cell_flag_index); 
*    tbox::Pointer<NodeData<DIM,int> > node_flag_data = 
*        patch->getPatchData(node_flag_index); 
*    tbox::Pointer<CellData<DIM,double> > vol_frac_data = 
*        patch->getPatchData(vol_frac_index); 
*    tbox::Pointer<IndexData<DIM,CutCell> > cut_cell_data = 
*        patch->getPatchData(cut_cell_index); 
*    

Input parameters specify the list of shapes to be used to construct the boundary, the desired accuracy of the volume and area fraction computation, and information about whether write constructed embedded boundary information to file, or read from file.

Required input keys and data types: NONE

Optional input keys, data types, and defaults:

-verbose boolean specifying whether to output information about the embedded boundary, such as the number of cut cells, the error in the volume calculation, etc. to the log file. If no value is supplied in input, the default is TRUE.

-max_subdivides
integer specifying the number of cell subdivides when computing the volume and area fractions. The larger the number of subdivides, the more accurate the fraction calculation will be, but the calculation will also be more expensive. If no value is supplied in input, a default of 0 is used.

-read_from_file bool specifying whether to read the embedded boundary from file.
If true, it will read flag, volume fraction, and cut cell information describing the embedded boundary from the specified HDF "dirname" (specified below). If no value is supplied, the default is FALSE.

-write_to_file bool specifying whether to write the embedded boundary to the specified HDF "dirname" (specified below). If true, it will write flag, volume fraction, and cut cell information computed while building the embedded boundary to the specified HDF "dirname" (specified below). If no value is supplied, the default is FALSE.

-dirname string specifying the name of the HDF directory to read/write embedded boundary information. Used with the "read_from_file" and "write_to_file" options.

-compute_areas_and_normal
bool specifying whether to compute areas and normal information.
Some applications only need the volume fraction so it is not necessary to invoke the extra cost of computing the areas and normal. If no value is supplied, the default is TRUE.

-compute_cutcell_index_data bool specifying whether to compute a list of cut cells. If false, it will compute the cell flag and volume fraction information but will not create and store the list of CutCell data structs. If no value is supplied, the default is TRUE.

-compute_boundary_node_data bool specifying whether to mark nodes that are just inside the geometry as BOUNDARY nodes, and to compute the centroid of the wetted cut area (i.e. the "front area") of the cut cells. This information may be needed for node-based finite difference or finite element computations of the embedded boundary. If no value is supplied, the default is TRUE.

-use_recursive_algs bool specifying whether to use a recursive algorithm to compute volume and area fractions. If true, it volume and area fractions are computed by recursively subdividing the cell until the max number of subdivides is reached. If false, it will apply the max subdivides to divide the cell into a (potentially large) array of subcells. Algorithmically, the recursive algorithm has fewer operations but the non-recursive algorithm may be more computationally efficient because it can do array-based operations. If no value is supplied, the default is TRUE.

-Shapes sub-database that specifies information about the analytic shapes used to construct an embedded boundary. See the EmbeddedBoundaryShapeSphere and EmbeddedBoundaryShapePolygon class headers for information on the inputs required.

-CubesPatchInterface sub-database for the Cubes interface, a cut-cell mesh generator from NASA Ames. See the CubesPatchInterface class header for information about the required inputs.

The following represents a sample input entry:

*  EmbeddedBoundaryGeometry{
*     max_subdivides = 2
*     read_from_file  = FALSE 
*     write_to_file   = FALSE
*     dirname         = "eb_grid"
*     compute_areas_and_normal   = TRUE 
*     compute_boundary_node_data = FALSE
*     use_recursive_algs         = FALSE 
*
*     Shapes {
*        Shape1 {
*           type = "POLYGON"    
*           vertices {
*              v1 = 1.0 , 1.0
*              v2 = 1.5 , .5 
*              v3 = 2.0 , 1.75 
*              v4 = 1.5 , 4.5
*              v5 = .5 , 2.0 
*           }
*           height = 5.0 // only used for 3D 
*        }
*        Shape2 {
*           type = "SPHERE"    
*           center = 65., 50.
*           radius = 20.     
*        }
*     }
*  }
*  

Note: Each shape has its own specific set of inputs. See the individual shapes for information about input requirements.

See also
appu::EmbeddedBoundaryShapeSphere
appu::EmbeddedBoundaryShapePolygon
appu::CubesPatchInterface
appu::ElevenPatchInterface
appu::CutCell

Member Enumeration Documentation

◆ CELL_TYPE

Enumerated type for the different cell classifications.

  • SOLID {Cell is located in the "solid" region.}
  • CUT {Cell is cut, meaning a CutCell data structure will be maintained at this cell.}
  • BORDER {Cell neighbors a cut cell, in the "flow" region.}
  • FLOW {Cell is located in the "flow" region.}
Enumerator
SOLID 
CUT 
BORDER 
FLOW 

◆ NODE_TYPE

Enumerated type for inside/outside node classification.

  • INSIDE {Node is located "inside" the prescribed geometry.}
  • OUTSIDE {Node is outside the geometry.}
  • BOUNDARY {Node is on the boundary of the geometry. That is it is the first one "inside" the geometry.}
  • ONBOUNDARY {Node is located exactly on the boundary of the geometry (used to avoid divide-by-zero problems) in numerical operations at embedded boundary.}
Enumerator
OUTSIDE 
INSIDE 
BOUNDARY 
ONBOUNDARY 

◆ PACK_RETURN_TYPE

  • ALL_ZERO (Fractions are 0.0 for every cell in patch.)
  • ALL_ONE (Fractions are 1.0 for every cell in patch.)
  • MIXTURE (There is some of this species/material in one or more cells of this patch, but the above two cases do not apply.)
Enumerator
VISIT_ALLZERO 
VISIT_ALLONE 
VISIT_MIXED 

Constructor & Destructor Documentation

◆ EmbeddedBoundaryGeometry()

template<int DIM>
SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::EmbeddedBoundaryGeometry ( const std::string &  object_name,
tbox::Pointer< tbox::Database input_db = NULL,
const tbox::Pointer< geom::CartesianGridGeometry< DIM > >  grid_geom = NULL,
const hier::IntVector< DIM > &  nghosts = hier::IntVector< DIM >(0) 
)

Constructor sets default values and reads data from input.

Parameters
object_nameName of object.
input_dbInput database.
grid_geomThe grid geometry (e.g. cartesian) used in the problem.
nghostsNumber of ghosts used to hold the embedded boundary. If not supplied, defaults to 0.

The grid_geom and nghosts may be NULL if the embedded boundary information is to be supplied by a file, either by restart or by other input such as CART3D. If the embedded boundary is to be constructed using analytic shapes in SAMRAI, the input_db and and grid_geom arguments must be supplied and may not be NULL.

◆ ~EmbeddedBoundaryGeometry()

Destructor deallocates data describing and unregisters the object with the restart manager if previously registered.

Member Function Documentation

◆ buildEmbeddedBoundaryOnLevel()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::buildEmbeddedBoundaryOnLevel ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy = NULL,
const tbox::Pointer< hier::PatchLevel< DIM > >  old_level = NULL 
)

Build an embedded boundary by forming the set of cut cells, extend them to appropriately apply the physical boundary conditions, and lastly compute the surrounding volumes for mass correction.

Depending on the arguments supplied, this method may be used in one of two ways. The first, taking only the level as an argument, performs an exhaustive search of all cells on the level to find the cut cells and classify the cells as inside or outside.
The second, taking as additional arguments a hierarchy, coarser_level, and possibly old_level, performs the same function but uses the information on the coarser and old levels to narrow the search for cut cells, making it considerably faster. Generally, the first method is used for the coarsest level only and the second is used for all subsequent finer levels.

Parameters
levelPatch level on which embedded boundary is to be constructed.
hierarchyPatch hierarchy of the level. Required if the level supplied in the first argument is in the hierarchy and is not the coarsest level.
old_levelPatch level which holds "old" embedded boundary data (not required). Use this if regridding and you have an embedded boundary at the old level to speed construction of the boundary on the new level.

◆ tagInsideOutsideNodesOnLevel()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::tagInsideOutsideNodesOnLevel ( const tbox::Pointer< hier::PatchLevel< DIM > >  level)

Tag nodes as being inside or outside the geometry. Some applications only require this knowledge and do not use need the volume/area fraction information for cut cells.

Parameters
levelPatch level where tagging takes place.

◆ registerVisItDataWriter()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::registerVisItDataWriter ( tbox::Pointer< appu::VisItDataWriter< DIM > >  visit_writer)

Register a VisIt data writer so this class will write plot files that may be postprocessed with the VisIt visualization tool.

Parameters
visit_writerVisIt data writer

◆ getCellFlagDataId()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getCellFlagDataId ( ) const

Return the descriptor index for the CellData<DIM,int> patch data that holds the integer flag at each cell, identifying it as being solid, cut, boundary, or flow.

◆ getCellVolumeDataId()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getCellVolumeDataId ( ) const

Return the descriptor index for the CellData<DIM,double> patch data that holds the cell volume fraction.

◆ getIndexCutCellDataId()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getIndexCutCellDataId ( ) const

Return the descriptor index for the IndexData<DIM,CutCell> patch data that holds the list of embedded boundary cells on the patch.

◆ getNodeInsideOutsideDataId()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getNodeInsideOutsideDataId ( ) const

Return the descriptor index for the NodeData<DIM,int> patch data that specifies a cell as being inside or outside the geometry.

◆ computeTotalVolumeOnLevel()

template<int DIM>
double SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::computeTotalVolumeOnLevel ( const tbox::Pointer< hier::PatchLevel< DIM > >  level)

Compute the total volume, which is the sum of the volume fractions, on the supplied level.

◆ setGridGeometry()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::setGridGeometry ( const tbox::Pointer< geom::CartesianGridGeometry< DIM > >  grid_geom)

Set the grid geometry, if it was not supplied when the object was constructed.

◆ writeLevelEmbeddedBoundaryDataToFile()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::writeLevelEmbeddedBoundaryDataToFile ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const std::string &  dirname 
) const

Write embedded boundary information - cell flag and cut cell information - to supplied directory. Files of the form "ebmesh-l<ln>-p<pid>.hdf", where ln is the level number and pid is the processor number, will be written to the directory.

◆ packMaterialFractionsIntoDoubleBuffer()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::packMaterialFractionsIntoDoubleBuffer ( double dbuffer,
const hier::Patch< DIM > &  patch,
const hier::Box< DIM > &  region,
const std::string &  material_name 
) const
virtual

The following routine:

packMaterialFractionsIntoDoubleBuffer()

is a concrete implementation of function declared in the appu::VisMaterialsDataStrategy abstract base class.

Put volume data located on the patch into the double buffer over the specified region.

Parameters
dbufferdouble buffer into which materials data is packed.
patchsupplied patch on which materials data is defined
regionregion over which data is packed
material_namename of the material

Reimplemented from SAMRAI::appu::VisMaterialsDataStrategy< DIM >.

◆ setPhysicalBoundaryConditions()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::setPhysicalBoundaryConditions ( hier::Patch< DIM > &  patch,
const double  fill_time,
const hier::IntVector< DIM > &  ghost_width_to_fill 
)
virtual

The following routines:

setPhysicalBoundaryConditions(),
getRefineOpStencilWidth(),
postprocessRefine()

are concrete implementations of functions declared in the xfer::RefinePatchStrategy abstract base class.

Set the data in ghost cells corresponding to physical boundary conditions. Specific boundary conditions are determined by information specified in input file and numerical routines.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ preprocessRefine()

template<int DIM>
virtual void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::preprocessRefine ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::Box< DIM > &  fine_box,
const hier::IntVector< DIM > &  ratio 
)
inlinevirtual

Perform user-defined refining operations. This member function is called before the other refining operators. For this class, no preprocessing is needed for the refine operators so it is setup to do nothing.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ postprocessRefine()

template<int DIM>
virtual void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::postprocessRefine ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::Box< DIM > &  fine_box,
const hier::IntVector< DIM > &  ratio 
)
virtual

Postprocess data after the refinement operator is applied. For the embedded boundary data, refining involves two steps; First, refine the "flag" data, which will designate on the fine level where the boundary exists; Second, on each fine cell that is flagged to possibly contain the boundary, go through and determine whether the cell is indeed cut by the boundary and adjust the volume, flag, and boundary cell data on the fine level as necessary. This method invokes step 2 of this refine operation.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ getRefineOpStencilWidth()

template<int DIM>
virtual hier::IntVector<DIM> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getRefineOpStencilWidth ( ) const
virtual

Return maximum stencil width needed for user-defined data interpolation operations. Default is to return zero, assuming no user-defined operations provided.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ putToDatabase()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::putToDatabase ( tbox::Pointer< tbox::Database db)
virtual

The following routine:

putToDatabase()

is a concrete implementation of a function declared in the tbox::Serializable abstract base class.

Implements SAMRAI::tbox::Serializable.

◆ initializeVariables()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::initializeVariables ( const hier::IntVector< DIM > &  nghosts)
private

◆ computeEmbeddedBoundaryOnLevelWithPackage()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::computeEmbeddedBoundaryOnLevelWithPackage ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
int cut_cells_on_level 
)
private

◆ computeEmbeddedBoundaryOnLevel()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::computeEmbeddedBoundaryOnLevel ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
int cut_cells_on_level,
double l2_volume_error,
double l2_area_error,
double max_volume_error,
double max_area_error 
)
private

◆ calculateCutCellInformation()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::calculateCutCellInformation ( appu::CutCell< DIM > &  cut_cell,
const double lower,
const double upper,
const double fullcellvol,
const double fullareas,
double volume_error_estimate,
double area_error_estimate 
)
private

◆ calculateBoundaryNodeInformation()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::calculateBoundaryNodeInformation ( appu::CutCell< DIM > &  cut_cell,
tbox::Pointer< pdat::NodeData< DIM, int > > &  node_flag,
const double lower,
const double upper,
const double dx 
)
private

◆ calculateVolume()

template<int DIM>
double SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::calculateVolume ( const double cell_lower,
const double cell_upper,
double error_estimate 
) const
private

◆ recursiveCalculateVolume()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::recursiveCalculateVolume ( const double cell_lower,
const double cell_upper,
double volume,
double error_estimate,
int subdivide_level 
) const
private

◆ calculateArea()

template<int DIM>
double SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::calculateArea ( const double cell_lower,
const double cell_upper,
const int  face,
double error_estimate 
) const
private

◆ recursiveCalculateArea()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::recursiveCalculateArea ( const double cell_lower,
const double cell_upper,
const int  face,
double area,
double error_estimate,
int subdivide_level 
) const
private

◆ classifyCell()

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::classifyCell ( const double lower,
const double upper 
) const
private

◆ refineEmbeddedBoundary()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::refineEmbeddedBoundary ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const tbox::Pointer< hier::PatchLevel< DIM > >  old_level = NULL,
const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy = NULL 
)
private

◆ setEmbeddedBoundaryAtPhysicalBoundaries()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::setEmbeddedBoundaryAtPhysicalBoundaries ( const tbox::Pointer< hier::PatchLevel< DIM > >  level)
private

Set the embedded boundary cells on physical boundarys according to the specified boundary conditions.

◆ setLevelRatioInformation()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::setLevelRatioInformation ( const tbox::Pointer< hier::PatchLevel< DIM > >  level)
private

◆ doNativeShapeInsideOutside()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::doNativeShapeInsideOutside ( const int nx,
const double dx,
const double origin,
int node_flag 
) const
private

◆ readLevelEmbeddedBoundaryDataFromFile()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::readLevelEmbeddedBoundaryDataFromFile ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const std::string &  dirname 
) const
private

Read embedded boundary information - cell flag and cut cell information - from supplied directory. Files of the form "ebmesh-l<ln>-p<pid>.hdf", where ln is the level number and pid is the processor number, should exist in the directory. If they don't exist, an error will be thrown.

◆ getFromInput()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getFromInput ( tbox::Pointer< tbox::Database db,
const bool  is_from_restart 
)
private

◆ getFromRestart()

template<int DIM>
void SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::getFromRestart ( )
private

◆ packMaterialFractionsIntoSparseBuffers()

template<int DIM>
virtual int SAMRAI::appu::VisMaterialsDataStrategy< DIM >::packMaterialFractionsIntoSparseBuffers ( int mat_list,
std::vector< int > &  mix_zones,
std::vector< int > &  mix_mat,
std::vector< double > &  vol_fracs,
std::vector< int > &  next_mat,
const hier::Patch< DIM > &  patch,
const hier::Box< DIM > &  region 
) const
virtualinherited

This function, which must be implemented whenever materials are used if the sparse packing format is to be used. Packing traverses the patch in a column major order (i.e. (f(x_0,y_0,z_0), f(x_1,y_0,z_0), f(x_2,y_0,z_0), ...) ). If the current cell is clean, then the material number of the material occupying the cell should be packed in to the buffer. If the cell is mixed, then a negative index (i.e., with numbering begining at -1) into the auxilliary mix_zones, mix_mat, vol_fracs, and next_mat vectors which will contain the sparse representation of the volume fractions (VisIt will correct the negative index internally and offset to a zero based representation). For each component of a mixed cell packMaterialFractionsIntoDoubleBuffer() should pack an entry in the following vectors: mix_zones: the cell number with which the fraction is associated mix_mat: the material number of the fraction vol_fracs: the volume fraction of the current material next_mat: either the (positive but still offset by one) index to the next mixed component within the auxilliary vectors or a 0, indicating the end of the mixed components for this cell.

If a non-zero ghost cell vector was specified when registerMaterialNames() was invoked, then ghost data corresponding to this ghost cell vector must be packed into this double buffer.

This method will be called once for each patch.

A enumerated PACK_RETURN_TYPE is used for a return value. To save space in the visit data file, you may choose to set the return value to ALL_ONE to indicate that a single material occupies 100% of each cell on the patch. Otherwise, MIXTURE should be returned even if individual cells are not mixed, as this indicates that the patch contains multiple materials.

Parameters
bufferDouble precision array into which the material numbers (or negative indices) are packed.
mix_zonesstd::vector<int> into which the cell number associated with the mixed components are packed.
mix_matstd::vector<int> into which the material numbers of the mixed components are packed.
vol_fracsstd::vector<double> into which the volume fractions (between 0.0 and 1.0) of the mixed components are packed.
next_matstd::vector<int> into which the (positive) index of the next mixed component or a terminating 0 are packed.
patchhier::Patch on which fractions are defined.
regionhier::Box region over which to pack fractions.
Returns
The enumeration constant VisMaterialsDataStrategy::ALL_ONE, or VisMaterialsDataStrategy::MIXTURE.

◆ packSpeciesFractionsIntoDoubleBuffer()

template<int DIM>
virtual int SAMRAI::appu::VisMaterialsDataStrategy< DIM >::packSpeciesFractionsIntoDoubleBuffer ( double buffer,
const hier::Patch< DIM > &  patch,
const hier::Box< DIM > &  region,
const std::string &  material_name,
const std::string &  species_name 
) const
virtualinherited

This user supplied function packs species fractions of the given material, patch, and region into the supplied 1D double precision buffer. If a non-zero ghost cell vector was specified when registerSpeciesNames() was invoked, then ghost data corresponding to this ghost cell vector must be packed into this double buffer. The data must be packed into the buffer in column major order.

This method will be called once for each species for each patch.

The method must return a PACK_RETURN_TYPE of ALL_ONE, ALL_ZERO, or MIXED. See the discussion above for the "packMaterialFractionsIntoDoubleBuffer()" method for an explanation of correct return values.

Parameters
bufferDouble precision array into which cell-centered species fractions are packed.
patchhier::Patch on which fractions are defined.
regionhier::Box region over which to pack fractions.
material_nameString identifier for the material to which the species belongs.
species_nameString identifier for the species.
Returns
The enumeration constant VisMaterialsDataStrategy::ALL_ZERO, VisMaterialsDataStrategy::ALL_ONE, or VisMaterialsDataStrategy::MIXTURE.

◆ preprocessRefineBoxes()

template<int DIM>
virtual void SAMRAI::xfer::RefinePatchStrategy< DIM >::preprocessRefineBoxes ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::BoxList< DIM > &  fine_boxes,
const hier::IntVector< DIM > &  ratio 
)
virtualinherited

Virtual function to perform user-defined refine operations. This member function is called before standard refining operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The preprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box regions.

Typically, only the pure virtual members of this class are implemented in user-defined subclasses of this base class. This version of the preprocess function operates on an entire box list. By default, this version simply loops over the box list and calls the single-box version, which is a pure virtual function.

Parameters
fineFine patch containing destination data.
coarseCoarse patch containing source data.
fine_boxestbox::List of box regions on fine patch into which data is refined.
ratioInteger vector containing ratio relating index space between coarse and fine patches.

Reimplemented in SAMRAI::solv::CartesianRobinBcHelper< DIM >.

◆ postprocessRefineBoxes()

template<int DIM>
virtual void SAMRAI::xfer::RefinePatchStrategy< DIM >::postprocessRefineBoxes ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::BoxList< DIM > &  fine_boxes,
const hier::IntVector< DIM > &  ratio 
)
virtualinherited

Virtual function to perform user-defined refine operations. This member function is called after standard refining operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The postprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box regions.

Typically, only the pure virtual members of this class are implemented in user-defined subclasses of this base class. This version of the postprocess function operates on an entire box list. By default, this version simply loops over the box list and calls the single-box version, which is a pure virtual function.

Parameters
fineFine patch containing destination data.
coarseCoarse patch containing source data.
fine_boxestbox::List of box regions on fine patch into which data is refined.
ratioInteger vector containing ratio relating index space between coarse and fine patches.

Reimplemented in SAMRAI::solv::CartesianRobinBcHelper< DIM >.

Member Data Documentation

◆ d_object_name

template<int DIM>
std::string SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_object_name
private

◆ d_grid_geometry

template<int DIM>
tbox::Pointer<geom::CartesianGridGeometry<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_grid_geometry
private

◆ d_visit_writer

template<int DIM>
tbox::Pointer<appu::VisItDataWriter<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_visit_writer
private

◆ d_verbose

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_verbose
private

◆ d_flow_type

template<int DIM>
std::string SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_flow_type
private

◆ d_ebdry_nghosts

template<int DIM>
hier::IntVector<DIM> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ebdry_nghosts
private

◆ d_ebdry_var

template<int DIM>
tbox::Pointer< pdat::IndexVariable<DIM,appu::CutCell<DIM>, pdat::CellGeometry<DIM> > > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ebdry_var
private

◆ d_cell_flag_var

template<int DIM>
tbox::Pointer< pdat::CellVariable<DIM,int> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cell_flag_var
private

◆ d_cell_vol_var

template<int DIM>
tbox::Pointer< pdat::CellVariable<DIM,double> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cell_vol_var
private

◆ d_node_flag_var

template<int DIM>
tbox::Pointer< pdat::NodeVariable<DIM,int> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_node_flag_var
private

◆ d_ebdry_data_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ebdry_data_id
private

◆ d_ebdry_scratch_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ebdry_scratch_id
private

◆ d_cell_flag_data_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cell_flag_data_id
private

◆ d_cell_flag_scratch_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cell_flag_scratch_id
private

◆ d_cell_vol_data_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cell_vol_data_id
private

◆ d_cell_vol_scratch_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cell_vol_scratch_id
private

◆ d_node_flag_data_id

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_node_flag_data_id
private

◆ d_ebdry_refine_alg

template<int DIM>
tbox::Pointer<xfer::RefineAlgorithm<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ebdry_refine_alg
private

◆ d_max_levels

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_max_levels
private

◆ d_ratio_to_coarser

template<int DIM>
tbox::Array<hier::IntVector<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ratio_to_coarser
private

◆ d_level_ratios

template<int DIM>
tbox::Array<hier::IntVector<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_level_ratios
private

◆ d_use_cubes

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_use_cubes
private

◆ d_cubes_interface

template<int DIM>
tbox::Pointer<appu::CubesPatchInterface<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_cubes_interface
private

◆ d_use_eleven_inside_outside

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_use_eleven_inside_outside
private

◆ d_use_eleven_boundary_node

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_use_eleven_boundary_node
private

◆ d_eleven_interface

template<int DIM>
tbox::Pointer<ElevenPatchInterface<DIM> > SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_eleven_interface
private

◆ d_shapes

◆ d_max_subdivides

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_max_subdivides
private

◆ d_read_ebdry

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_read_ebdry
private

◆ d_write_ebdry

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_write_ebdry
private

◆ d_ebdry_dirname

template<int DIM>
std::string SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_ebdry_dirname
private

◆ d_use_recursive_algs

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_use_recursive_algs
private

◆ d_compute_areas_and_normal

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_compute_areas_and_normal
private

◆ d_compute_cutcell_index_data

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_compute_cutcell_index_data
private

◆ d_compute_boundary_node_data

template<int DIM>
bool SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_compute_boundary_node_data
private

◆ d_node_bdry_cond

template<int DIM>
tbox::Array<int> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_node_bdry_cond
private

◆ d_edge_bdry_cond

template<int DIM>
tbox::Array<int> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_edge_bdry_cond
private

◆ d_face_bdry_cond

template<int DIM>
tbox::Array<int> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_face_bdry_cond
private

◆ t_compute_eb

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_compute_eb
private

◆ t_calc_node_inout

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_calc_node_inout
private

◆ t_calc_volume

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_calc_volume
private

◆ t_calc_area

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_calc_area
private

◆ t_calc_boundary_node

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_calc_boundary_node
private

◆ t_eb_phys_bdry

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_eb_phys_bdry
private

◆ t_read_write_eb

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::t_read_write_eb
private

◆ d_step_count

template<int DIM>
int SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::d_step_count
private

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