SAMRAI::appu::EmbeddedBoundaryGeometry< DIM > Class Template Reference

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

#include <source/apputils/embedded_boundary/EmbeddedBoundaryGeometry.h>

Inheritance diagram for SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >:

Inheritance graph
[legend]
List of all members.

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 }

Public Member Functions

 EmbeddedBoundaryGeometry (const std::string &object_name, tbox::Pointer< tbox::Database > input_db=(0), const tbox::Pointer< geom::CartesianGridGeometry< DIM > > grid_geom=(0), 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=(0), const tbox::Pointer< hier::PatchLevel< DIM > > old_level=(0))
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)

Detailed Description

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

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

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.

2. 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.

3. Access information about the embedded boundary from patches:

 *    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

template<int DIM>
enum SAMRAI::appu::EmbeddedBoundaryGeometry::CELL_TYPE

Enumerated type for the different cell classifications.

Enumerator:
SOLID 
CUT 
BORDER 
FLOW 

template<int DIM>
enum SAMRAI::appu::EmbeddedBoundaryGeometry::NODE_TYPE

Enumerated type for inside/outside node classification.

Enumerator:
OUTSIDE 
INSIDE 
BOUNDARY 
ONBOUNDARY 


Constructor & Destructor Documentation

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

Constructor sets default values and reads data from input.

Parameters:
object_name Name of object.
input_db Input database.
grid_geom The grid geometry (e.g. cartesian) used in the problem.
nghosts Number 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.

template<int DIM>
SAMRAI::appu::EmbeddedBoundaryGeometry< DIM >::~EmbeddedBoundaryGeometry< DIM > (  ) 

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


Member Function Documentation

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

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:
level Patch level on which embedded boundary is to be constructed.
hierarchy Patch hierarchy of the level. Required if the level supplied in the first argument is in the hierarchy and is not the coarsest level.
old_level Patch 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.

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:
level Patch level where tagging takes place.

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_writer VisIt data writer

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.

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.

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.

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.

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.

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.

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.

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]

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

Parameters:
dbuffer double buffer into which materials data is packed.
patch supplied patch on which materials data is defined
region region over which data is packed
material_name name of the material

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

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]

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 >.

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 
) [inline, virtual]

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 >.

template<int DIM>
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 >.

template<int DIM>
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 >.

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.


The documentation for this class was generated from the following files:
Generated on Thu Jun 18 11:28:21 2009 for SAMRAI by  doxygen 1.5.1