DOLFIN
DOLFIN C++ interface
|
#include <MultiMesh.h>
Public Types | |
typedef std::pair< std::vector< double >, std::vector< double > > | quadrature_rule |
Structure storing a quadrature rule. | |
typedef std::vector< Point > | Simplex |
A simplex is a list of points. | |
typedef std::pair< std::vector< Simplex >, std::set< std::size_t > > | Polyhedron |
A polyhedron is a list of simplices and the part numbers. | |
typedef std::vector< std::size_t > | IncExcKey |
Key to identify polyhedra. | |
Public Member Functions | |
MultiMesh () | |
Create empty multimesh. | |
MultiMesh (std::vector< std::shared_ptr< const Mesh >> meshes, std::size_t quadrature_order) | |
Create multimesh from given list of meshes. | |
MultiMesh (std::shared_ptr< const Mesh > mesh_0, std::size_t quadrature_order) | |
Create multimesh from one mesh. | |
MultiMesh (std::shared_ptr< const Mesh > mesh_0, std::shared_ptr< const Mesh > mesh_1, std::size_t quadrature_order) | |
Create multimesh from two meshes. | |
MultiMesh (std::shared_ptr< const Mesh > mesh_0, std::shared_ptr< const Mesh > mesh_1, std::shared_ptr< const Mesh > mesh_2, std::size_t quadrature_order) | |
Create multimesh from three meshes. | |
~MultiMesh () | |
Destructor. | |
std::size_t | num_parts () const |
std::shared_ptr< const Mesh > | part (std::size_t i) const |
const std::vector< unsigned int > & | uncut_cells (std::size_t part) const |
const std::vector< unsigned int > | cut_cells (std::size_t part) const |
const std::vector< unsigned int > & | covered_cells (std::size_t part) const |
void | mark_covered (std::size_t part, const std::vector< unsigned int > &cells) |
const std::map< unsigned int, std::vector< std::pair< std::size_t, unsigned int > > > & | collision_map_cut_cells (std::size_t part) const |
const std::map< unsigned int, quadrature_rule > & | quadrature_rules_cut_cells (std::size_t part) const |
const quadrature_rule | quadrature_rules_cut_cells (std::size_t part, unsigned int cell_index) const |
const std::map< unsigned int, std::vector< quadrature_rule > > & | quadrature_rules_overlap (std::size_t part) const |
const std::vector< quadrature_rule > | quadrature_rules_overlap (std::size_t part, unsigned int cell) const |
const std::map< unsigned int, std::vector< quadrature_rule > > & | quadrature_rules_interface (std::size_t part) const |
const std::vector< quadrature_rule > | quadrature_rules_interface (std::size_t part, unsigned int cell_index) const |
const std::map< unsigned int, std::vector< std::vector< double > > > & | facet_normals (std::size_t part) const |
std::shared_ptr< const BoundingBoxTree > | bounding_box_tree (std::size_t part) const |
std::shared_ptr< const BoundingBoxTree > | bounding_box_tree_boundary (std::size_t part) const |
void | add (std::shared_ptr< const Mesh > mesh) |
void | build (std::size_t quadrature_order=2) |
Build multimesh. | |
bool | is_built () const |
Check whether multimesh has been built. | |
void | clear () |
Clear multimesh. | |
double | compute_area () const |
double | compute_volume () const |
Corresponding function for volume. | |
std::string | plot_matplotlib (double delta_z=1, const std::string &filename="") const |
Create matplotlib string to plot 2D multimesh (small meshes only) | |
void | auto_cover (std::size_t p, const Point &point) |
Public Member Functions inherited from dolfin::Variable | |
Variable () | |
Create unnamed variable. | |
Variable (const std::string name, const std::string label) | |
Create variable with given name and label. | |
Variable (const Variable &variable) | |
Copy constructor. | |
virtual | ~Variable () |
Destructor. | |
const Variable & | operator= (const Variable &variable) |
Assignment operator. | |
void | rename (const std::string name, const std::string label) |
Rename variable. | |
std::string | name () const |
Return name. | |
std::string | label () const |
Return label (description) | |
std::size_t | id () const |
virtual std::string | str (bool verbose) const |
Return informal string representation (pretty-print) | |
Static Public Member Functions | |
static Parameters | default_parameters () |
Default parameter values. | |
Additional Inherited Members | |
Public Attributes inherited from dolfin::Variable | |
Parameters | parameters |
Parameters. | |
This class represents a collection of meshes with arbitrary overlaps. A multimesh may be created from a set of standard meshes spaces by repeatedly calling add(), followed by a call to build(). Note that a multimesh is not useful until build() has been called.
void MultiMesh::auto_cover | ( | std::size_t | p, |
const Point & | point | ||
) |
Marks all uncut and cut cells connected to the given point as covered. This can be used for instance to mark a hole as covered where one point inside the hole is known.
std::shared_ptr< const BoundingBoxTree > MultiMesh::bounding_box_tree | ( | std::size_t | part | ) | const |
Return the bounding box tree for the mesh of the given part
Arguments part (std::size_t) The part number
Returns std::shared_ptr<const BoundingBoxTree> The bounding box tree
std::shared_ptr< const BoundingBoxTree > MultiMesh::bounding_box_tree_boundary | ( | std::size_t | part | ) | const |
Return the bounding box tree for the boundary mesh of the given part
Arguments part (std::size_t) The part number
Returns std::shared_ptr<const BoundingBoxTree> The bounding box tree
const std::map< unsigned int, std::vector< std::pair< std::size_t, unsigned int > > > & MultiMesh::collision_map_cut_cells | ( | std::size_t | part | ) | const |
Return the collision map for cut cells of the given part
Arguments part (std::size_t) The part number
Returns std::map<unsigned int, std::vector<std::pair<std::size_t, unsigned int> > > A map from cell indices of cut cells to a list of cutting cells. Each cutting cell is represented as a pair (part_number, cutting_cell_index).
double MultiMesh::compute_area | ( | ) | const |
Compute total interface area or the total volume of multimesh by summing up quadrature weights. If the area or volume of the domain mesh is known, this is a good test to verify that the mesh-mesh intersections and quadrature are correct.
const std::vector< unsigned int > & MultiMesh::covered_cells | ( | std::size_t | part | ) | const |
Return the list of covered cells for given part. The covered cells are defined as all cells that collide with the domain of any part with higher part number, but not with the boundary of that part; in other words cells that are completely covered by any other part (and which therefore are inactive).
Arguments part (std::size_t) The part number
Returns std::vector<unsigned int> List of covered cell indices for given part
const std::vector< unsigned int > MultiMesh::cut_cells | ( | std::size_t | part | ) | const |
Return the list of cut cells for given part. The cut cells are defined as all cells that collide with the boundary of any part with higher part number.
FIXME: Figure out whether this makes sense; a cell may collide with the boundary of part j but may still be covered completely by the domain of part j + 1. Possible solution is to for each part i check overlapping parts starting from the top and working back down to i + 1.
Arguments part (std::size_t) The part number
Returns std::vector<unsigned int> List of cut cell indices for given part
const std::map< unsigned int, std::vector< std::vector< double > > > & MultiMesh::facet_normals | ( | std::size_t | part | ) | const |
Return facet normals for the interface on the given part
Arguments part (std::size_t) The part number
Returns std::map<unsigned int, std::vector<std::vector<double> > > A map from cell indices of cut cells to facet normals on an interface part cutting through the cell. A separate list of facet normals, one for each quadrature point, is given for each cutting cell and stored in the same order as in the collision map. The facet normals for each set of quadrature points is stored as a contiguous flattened array, the length of which should be equal to the number of quadrature points multiplied by the geometric dimension. Puh!
void MultiMesh::mark_covered | ( | std::size_t | part, |
const std::vector< unsigned int > & | cells | ||
) |
Mark a set of cells as covered in the mesh.
Arguments part (std::size_t) The part number cells (std::vector<unsigned int>) The cells to be covered
std::size_t MultiMesh::num_parts | ( | ) | const |
Return the number of meshes (parts) of the multimesh
Returns std::size_t The number of meshes (parts) of the multimesh.
std::shared_ptr< const Mesh > MultiMesh::part | ( | std::size_t | i | ) | const |
const std::map< unsigned int, MultiMesh::quadrature_rule > & MultiMesh::quadrature_rules_cut_cells | ( | std::size_t | part | ) | const |
Return quadrature rules for cut cells on the given part
Arguments part (std::size_t) The part number
Returns std::map<unsigned int, std::pair<std::vector<double>, std::vector<double> > > A map from cell indices of cut cells to quadrature rules. Each quadrature rule is represented as a pair of a flattened array of quadrature points and a corresponding array of quadrature weights.
const MultiMesh::quadrature_rule MultiMesh::quadrature_rules_cut_cells | ( | std::size_t | part, |
unsigned int | cell_index | ||
) | const |
Return quadrature rule for a given cut cell on the given part
Arguments part (std::size_t) The part number cell (unsigned int) The cell index
Returns std::pair<std::vector<double>, std::vector<double> > A quadrature rule represented as a pair of a flattened array of quadrature points and a corresponding array of quadrature weights. An error is raised if the given cell is not in the map.
const std::map< unsigned int, std::vector< MultiMesh::quadrature_rule > > & MultiMesh::quadrature_rules_interface | ( | std::size_t | part | ) | const |
Return quadrature rules for the interface on the given part
Arguments part (std::size_t) The part number
Returns std::map<unsigned int, std::vector<std::pair<std::vector<double>, std::vector<double> > > > A map from cell indices of cut cells to quadrature rules on an interface part cutting through the cell. A separate quadrature rule is given for each cutting cell and stored in the same order as in the collision map. Each quadrature rule is represented as a pair of an array of quadrature points and a corresponding flattened array of quadrature weights.
const std::vector< MultiMesh::quadrature_rule > MultiMesh::quadrature_rules_interface | ( | std::size_t | part, |
unsigned int | cell_index | ||
) | const |
Return quadrature rules for the interface of a given cut cell on the given part
Arguments part (std::size_t) The part number cell (unsigned int) The cell index
Returns std::vector<std::pair<std::vector<double>, std::vector<double> > > A vector of quadrature rules on the cut cell. A separate quadrature rule is given for each cutting cell and stored in the same order as in the collision map. Each quadrature rule represented as a pair of a flattened array of quadrature points and a corresponding array of quadrature weights. An error is raised if the given cell is not in the map.
Developer note: this function is mainly useful from Python and could be replaced by a suitable typemap that would make the previous more general function accessible from Python.
const std::map< unsigned int, std::vector< MultiMesh::quadrature_rule > > & MultiMesh::quadrature_rules_overlap | ( | std::size_t | part | ) | const |
Return quadrature rules for the overlap on the given part.
Arguments part (std::size_t) The part number
Returns std::map<unsigned int, std::vector<std::pair<std::vector<double>, std::vector<double> > > > A map from cell indices of cut cells to quadrature rules. A separate quadrature rule is given for each cutting cell and stored in the same order as in the collision map. Each quadrature rule is represented as a pair of an array of quadrature points and a corresponding flattened array of quadrature weights.
const std::vector< MultiMesh::quadrature_rule > MultiMesh::quadrature_rules_overlap | ( | std::size_t | part, |
unsigned int | cell | ||
) | const |
Return quadrature rules for the overlap for a given cell on the given part.
Arguments part (std::size_t) The part number Returns std::vector<std::pair<std::vector<double>, std::vector<double> > > A vector of quadrature rules on the cut cell. A separate quadrature rule is given for each cutting cell and stored in the same order as in the collision map. A quadrature rule represented as a pair of a flattened array of quadrature points and a corresponding array of quadrature weights. An error is raised if the given cell is not in the map.
const std::vector< unsigned int > & MultiMesh::uncut_cells | ( | std::size_t | part | ) | const |
Return the list of uncut cells for given part. The uncut cells are defined as all cells that don't collide with any cells in any other part with higher part number.
Arguments part (std::size_t) The part number
Returns std::vector<unsigned int> List of uncut cell indices for given part