oomph::TetMeshFacetedSurface Class Reference

#include <tet_mesh.h>

+ Inheritance diagram for oomph::TetMeshFacetedSurface:

Public Member Functions

 TetMeshFacetedSurface ()
 Constructor: More...
 
virtual ~TetMeshFacetedSurface ()
 Empty destructor. More...
 
unsigned nvertex () const
 Number of vertices. More...
 
unsigned nfacet () const
 Number of facets. More...
 
unsigned one_based_facet_boundary_id (const unsigned &j) const
 One-based boundary id of j-th facet. More...
 
unsigned one_based_vertex_boundary_id (const unsigned &j) const
 First (of possibly multiple) one-based boundary id of j-th vertex. More...
 
double vertex_coordinate (const unsigned &j, const unsigned &i) const
 i-th coordinate of j-th vertex More...
 
unsigned nvertex_on_facet (const unsigned &j) const
 Number of vertices defining the j-th facet. More...
 
bool boundaries_can_be_split_in_tetgen ()
 Test whether boundary can be split in tetgen. More...
 
void enable_boundaries_can_be_split_in_tetgen ()
 Test whether boundaries can be split in tetgen. More...
 
void disable_boundaries_can_be_split_in_tetgen ()
 Test whether boundaries can be split in tetgen. More...
 
TetMeshFacetfacet_pt (const unsigned &j) const
 Pointer to j-th facet. More...
 
TetMeshVertexvertex_pt (const unsigned &j) const
 Pointer to j-th vertex. More...
 
DiskLikeGeomObjectWithBoundariesgeom_object_with_boundaries_pt ()
 
void output (std::ostream &outfile) const
 Output. More...
 
void output (const std::string &filename) const
 Output. More...
 
virtual void boundary_zeta01 (const unsigned &facet_id, const double &zeta_boundary, Vector< double > &zeta)
 
virtual void boundary_zeta12 (const unsigned &facet_id, const double &zeta_boundary, Vector< double > &zeta)
 
virtual void boundary_zeta20 (const unsigned &facet_id, const double &zeta_boundary, Vector< double > &zeta)
 
Vector< unsignedvertex_index_in_tetgen (const unsigned &f)
 

Protected Attributes

Vector< TetMeshVertex * > Vertex_pt
 Vector pointers to vertices. More...
 
Vector< TetMeshFacet * > Facet_pt
 Vector of pointers to facets. More...
 
bool Boundaries_can_be_split_in_tetgen
 
Vector< Vector< unsigned > > Facet_vertex_index_in_tetgen
 
DiskLikeGeomObjectWithBoundariesGeom_object_with_boundaries_pt
 GeomObject with boundaries associated with this surface. More...
 

Private Member Functions

void setup_facet_connectivity_for_tetgen ()
 Setup facet connectivity for tetgen. More...
 

Detailed Description

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// Base class for tet mesh boundary defined by polygonal planar facets

Constructor & Destructor Documentation

◆ TetMeshFacetedSurface()

oomph::TetMeshFacetedSurface::TetMeshFacetedSurface ( )
inline

Constructor:

312  {
313  }
bool Boundaries_can_be_split_in_tetgen
Definition: tet_mesh.h:490
DiskLikeGeomObjectWithBoundaries * Geom_object_with_boundaries_pt
GeomObject with boundaries associated with this surface.
Definition: tet_mesh.h:497

◆ ~TetMeshFacetedSurface()

virtual oomph::TetMeshFacetedSurface::~TetMeshFacetedSurface ( )
inlinevirtual

Empty destructor.

316 {}

Member Function Documentation

◆ boundaries_can_be_split_in_tetgen()

bool oomph::TetMeshFacetedSurface::boundaries_can_be_split_in_tetgen ( )
inline

Test whether boundary can be split in tetgen.

356  {
358  }

References Boundaries_can_be_split_in_tetgen.

Referenced by oomph::TetgenMesh< ELEMENT >::TetgenMesh().

◆ boundary_zeta01()

virtual void oomph::TetMeshFacetedSurface::boundary_zeta01 ( const unsigned facet_id,
const double zeta_boundary,
Vector< double > &  zeta 
)
inlinevirtual

Virtual function that specifies the variation of the zeta coordinates in the GeomObject along the boundary connecting vertices 0 and 1 in facet facet_id. Default implementation: Linear interpolation between edges. zeta_boundary=0.0: we're on vertex 0; zeta_boundary=1.0: we're on vertex 1.

Reimplemented in oomph::DiskTetMeshFacetedSurface.

420  {
421  Vector<Vector<double>> zeta_vertex(2);
422  zeta_vertex[0] = Facet_pt[facet_id]->vertex_pt(0)->zeta_in_geom_object();
423  zeta_vertex[1] = Facet_pt[facet_id]->vertex_pt(1)->zeta_in_geom_object();
424  zeta[0] = zeta_vertex[0][0] +
425  (zeta_vertex[1][0] - zeta_vertex[0][0]) * zeta_boundary;
426  zeta[1] = zeta_vertex[0][1] +
427  (zeta_vertex[1][1] - zeta_vertex[0][1]) * zeta_boundary;
428  }
Vector< TetMeshFacet * > Facet_pt
Vector of pointers to facets.
Definition: tet_mesh.h:486
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152

References Facet_pt, and Eigen::zeta().

Referenced by oomph::TetMeshBase::snap_nodes_onto_geometric_objects().

◆ boundary_zeta12()

virtual void oomph::TetMeshFacetedSurface::boundary_zeta12 ( const unsigned facet_id,
const double zeta_boundary,
Vector< double > &  zeta 
)
inlinevirtual

Virtual function that specifies the variation of the zeta coordinates in the GeomObject along the boundary connecting vertices 1 and 2 in facet facet_id. Default implementation: Linear interpolation between edges. zeta_boundary=0.0: we're on vertex 1; zeta_boundary=1.0: we're on vertex 2.

439  {
440  Vector<Vector<double>> zeta_vertex(2);
441  zeta_vertex[0] = Facet_pt[facet_id]->vertex_pt(1)->zeta_in_geom_object();
442  zeta_vertex[1] = Facet_pt[facet_id]->vertex_pt(2)->zeta_in_geom_object();
443  zeta[0] = zeta_vertex[0][0] +
444  (zeta_vertex[1][0] - zeta_vertex[0][0]) * zeta_boundary;
445  zeta[1] = zeta_vertex[0][1] +
446  (zeta_vertex[1][1] - zeta_vertex[0][1]) * zeta_boundary;
447  }

References Facet_pt, and Eigen::zeta().

Referenced by oomph::TetMeshBase::snap_nodes_onto_geometric_objects().

◆ boundary_zeta20()

virtual void oomph::TetMeshFacetedSurface::boundary_zeta20 ( const unsigned facet_id,
const double zeta_boundary,
Vector< double > &  zeta 
)
inlinevirtual

Virtual function that specifies the variation of the zeta coordinates in the GeomObject along the boundary connecting vertices 2 and 0 in facet facet_id. Default implementation: Linear interpolation between edges. zeta_boundary=0.0: we're on vertex 2; zeta_boundary=1.0: we're on vertex 0.

458  {
459  Vector<Vector<double>> zeta_vertex(2);
460  zeta_vertex[0] = Facet_pt[facet_id]->vertex_pt(2)->zeta_in_geom_object();
461  zeta_vertex[1] = Facet_pt[facet_id]->vertex_pt(0)->zeta_in_geom_object();
462  zeta[0] = zeta_vertex[0][0] +
463  (zeta_vertex[1][0] - zeta_vertex[0][0]) * zeta_boundary;
464  zeta[1] = zeta_vertex[0][1] +
465  (zeta_vertex[1][1] - zeta_vertex[0][1]) * zeta_boundary;
466  }

References Facet_pt, and Eigen::zeta().

Referenced by oomph::TetMeshBase::snap_nodes_onto_geometric_objects().

◆ disable_boundaries_can_be_split_in_tetgen()

void oomph::TetMeshFacetedSurface::disable_boundaries_can_be_split_in_tetgen ( )
inline

Test whether boundaries can be split in tetgen.

368  {
370  }

References Boundaries_can_be_split_in_tetgen.

◆ enable_boundaries_can_be_split_in_tetgen()

void oomph::TetMeshFacetedSurface::enable_boundaries_can_be_split_in_tetgen ( )
inline

Test whether boundaries can be split in tetgen.

362  {
364  }

References Boundaries_can_be_split_in_tetgen.

◆ facet_pt()

TetMeshFacet* oomph::TetMeshFacetedSurface::facet_pt ( const unsigned j) const
inline

Pointer to j-th facet.

374  {
375  return Facet_pt[j];
376  }
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Facet_pt, and j.

Referenced by oomph::GmshTetMesh< ELEMENT >::build_it(), oomph::TetgenMesh< ELEMENT >::setup_reverse_lookup_schemes_for_faceted_surface(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), and oomph::GmshTetScaffoldMesh::write_geo_file().

◆ geom_object_with_boundaries_pt()

DiskLikeGeomObjectWithBoundaries* oomph::TetMeshFacetedSurface::geom_object_with_boundaries_pt ( )
inline

Access to GeomObject with boundaries associated with this surface (Null if there isn't one!)

387  {
389  }

References Geom_object_with_boundaries_pt.

Referenced by oomph::TetMeshBase::snap_nodes_onto_geometric_objects().

◆ nfacet()

◆ nvertex()

unsigned oomph::TetMeshFacetedSurface::nvertex ( ) const
inline

Number of vertices.

320  {
321  return Vertex_pt.size();
322  }
Vector< TetMeshVertex * > Vertex_pt
Vector pointers to vertices.
Definition: tet_mesh.h:483

References Vertex_pt.

Referenced by oomph::TetgenMesh< ELEMENT >::build_tetgenio(), oomph::GmshTetScaffoldMesh::write_geo_file(), and oomph::TetMeshFacetedClosedSurfaceForRemesh::~TetMeshFacetedClosedSurfaceForRemesh().

◆ nvertex_on_facet()

unsigned oomph::TetMeshFacetedSurface::nvertex_on_facet ( const unsigned j) const
inline

Number of vertices defining the j-th facet.

350  {
351  return Facet_pt[j]->nvertex();
352  }

References Facet_pt, and j.

◆ one_based_facet_boundary_id()

◆ one_based_vertex_boundary_id()

unsigned oomph::TetMeshFacetedSurface::one_based_vertex_boundary_id ( const unsigned j) const
inline

First (of possibly multiple) one-based boundary id of j-th vertex.

338  {
339  return Vertex_pt[j]->one_based_boundary_id();
340  }

References j, and Vertex_pt.

Referenced by oomph::TetgenMesh< ELEMENT >::build_tetgenio().

◆ output() [1/2]

void oomph::TetMeshFacetedSurface::output ( const std::string &  filename) const
inline

Output.

404  {
405  std::ofstream outfile;
406  outfile.open(filename.c_str());
407  output(outfile);
408  outfile.close();
409  }
void output(std::ostream &outfile) const
Output.
Definition: tet_mesh.h:392
string filename
Definition: MergeRestartFiles.py:39

References MergeRestartFiles::filename, and output().

◆ output() [2/2]

void oomph::TetMeshFacetedSurface::output ( std::ostream &  outfile) const
inline

Output.

393  {
394  unsigned n = Facet_pt.size();
395  for (unsigned j = 0; j < n; j++)
396  {
397  Facet_pt[j]->output(outfile);
398  }
399  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11

References Facet_pt, j, and n.

Referenced by FallingBlockProblem< ELEMENT >::FallingBlockProblem(), output(), and TetmeshPoissonProblem< ELEMENT >::TetmeshPoissonProblem().

◆ setup_facet_connectivity_for_tetgen()

void oomph::TetMeshFacetedSurface::setup_facet_connectivity_for_tetgen ( )
inlineprivate

Setup facet connectivity for tetgen.

503  {
504  unsigned nv_overall = Vertex_pt.size();
505  unsigned nf = nfacet();
506  Facet_vertex_index_in_tetgen.resize(nf);
507  for (unsigned f = 0; f < nf; f++)
508  {
509  unsigned nv = Facet_pt[f]->nvertex();
510  Facet_vertex_index_in_tetgen[f].resize(nv);
511  for (unsigned v = 0; v < nv; v++)
512  {
513  TetMeshVertex* my_vertex_pt = Facet_pt[f]->vertex_pt(v);
514  for (unsigned j = 0; j < nv_overall; j++)
515  {
516  if (my_vertex_pt == Vertex_pt[j])
517  {
519  break;
520  }
521  }
522  }
523  }
524  }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
unsigned nfacet() const
Number of facets.
Definition: tet_mesh.h:325
Vector< Vector< unsigned > > Facet_vertex_index_in_tetgen
Definition: tet_mesh.h:494
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237

References f(), Facet_pt, Facet_vertex_index_in_tetgen, j, nfacet(), v, and Vertex_pt.

Referenced by vertex_index_in_tetgen().

◆ vertex_coordinate()

double oomph::TetMeshFacetedSurface::vertex_coordinate ( const unsigned j,
const unsigned i 
) const
inline

i-th coordinate of j-th vertex

344  {
345  return Vertex_pt[j]->x(i);
346  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i, j, and Vertex_pt.

Referenced by oomph::TetgenMesh< ELEMENT >::build_tetgenio().

◆ vertex_index_in_tetgen()

Vector<unsigned> oomph::TetMeshFacetedSurface::vertex_index_in_tetgen ( const unsigned f)
inline

Facet connectivity: vertex_index[j] is the index of the j-th vertex (in the Vertex_pt vector) in facet f. Bit of an obscure functionality that's only needed for setup tetgen_io.

473  {
474  if (Facet_vertex_index_in_tetgen.size() != nfacet())
475  {
477  }
479  }
void setup_facet_connectivity_for_tetgen()
Setup facet connectivity for tetgen.
Definition: tet_mesh.h:502

References f(), Facet_vertex_index_in_tetgen, nfacet(), and setup_facet_connectivity_for_tetgen().

Referenced by oomph::TetgenMesh< ELEMENT >::build_tetgenio().

◆ vertex_pt()

TetMeshVertex* oomph::TetMeshFacetedSurface::vertex_pt ( const unsigned j) const
inline

Pointer to j-th vertex.

380  {
381  return Vertex_pt[j];
382  }

References j, and Vertex_pt.

Referenced by oomph::DiskTetMeshFacetedSurface::boundary_zeta01(), and oomph::GmshTetScaffoldMesh::write_geo_file().

Member Data Documentation

◆ Boundaries_can_be_split_in_tetgen

bool oomph::TetMeshFacetedSurface::Boundaries_can_be_split_in_tetgen
protected

Boolean to indicate whether extra vertices can be added on the boundary in tetgen

Referenced by boundaries_can_be_split_in_tetgen(), disable_boundaries_can_be_split_in_tetgen(), and enable_boundaries_can_be_split_in_tetgen().

◆ Facet_pt

◆ Facet_vertex_index_in_tetgen

Vector<Vector<unsigned> > oomph::TetMeshFacetedSurface::Facet_vertex_index_in_tetgen
protected

Facet connectivity: Facet_vertex_index[f][j] is the index of the j-th vertex (in the Vertex_pt vector) in facet f.

Referenced by setup_facet_connectivity_for_tetgen(), and vertex_index_in_tetgen().

◆ Geom_object_with_boundaries_pt

◆ Vertex_pt


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