oomph::MeshAsGeomObject Class Reference

#include <mesh_as_geometric_object.h>

+ Inheritance diagram for oomph::MeshAsGeomObject:

Public Member Functions

SamplePointContainersample_point_container_pt () const
 Pointer to the sample point container. More...
 
FiniteElementfinite_element_pt (const unsigned &e)
 Return pointer to e-th finite element. More...
 
unsigned sample_point_container_version () const
 
unsigned nelement ()
 Number of elements in the underlying mesh. More...
 
 MeshAsGeomObject (Mesh *const &mesh_pt)
 Constructor. More...
 
 MeshAsGeomObject (SamplePointContainerParameters *sample_point_container_parameters_pt)
 Constructor. More...
 
 MeshAsGeomObject ()
 Empty Constructor. More...
 
 ~MeshAsGeomObject ()
 Destructor. More...
 
 MeshAsGeomObject (const MeshAsGeomObject &)=delete
 Broken copy constructor. More...
 
void operator= (const MeshAsGeomObject &)=delete
 Broken assignment operator. More...
 
unsigned ngeom_data () const
 How many items of Data does the shape of the object depend on? More...
 
Datageom_data_pt (const unsigned &j)
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition (const Vector< double > &xi, DenseMatrix< double > &drdxi) const
 Return the derivative of the position. More...
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition_dt (const Vector< double > &zeta, const unsigned &j, Vector< double > &drdt)
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 

Private Member Functions

void build_it (SamplePointContainerParameters *sample_point_container_parameters_pt)
 Helper function to actually build the thing. More...
 

Private Attributes

Vector< Data * > Geom_data_pt
 Vector of pointers to Data items that affects the object's shape. More...
 
Vector< FiniteElement * > Sub_geom_object_pt
 Internal storage for the elements that constitute the object. More...
 
SamplePointContainerSample_point_container_pt
 Pointer to the sample point container. More...
 
MeshMesh_pt
 Pointer to mesh. More...
 
unsigned Sample_point_container_version
 

Additional Inherited Members

- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Detailed Description

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// This class provides a GeomObject representation of a given finite element mesh. The Lagrangian coordinate is taken to be the dimension of the (first) element in the mesh and the Eulerian coordinate is taken to be the dimension of the (first) node in the mesh. If there are no elements or nodes the appropriate dimensions will be set to zero. The constituent elements of the mesh must have their own GeomObject representations, so they must be FiniteElements, and they become sub-objects in this compound GeomObject.

Constructor & Destructor Documentation

◆ MeshAsGeomObject() [1/4]

oomph::MeshAsGeomObject::MeshAsGeomObject ( Mesh *const &  mesh_pt)
inline

Constructor.

314  : GeomObject()
315  {
316  // Create default parameters
317  SamplePointContainerParameters* sample_point_container_parameters_pt = 0;
319  mesh_pt, sample_point_container_parameters_pt);
320 
321  // Build the bastard
322  build_it(sample_point_container_parameters_pt);
323  delete sample_point_container_parameters_pt;
324  }
GeomObject()
Default constructor.
Definition: geom_objects.h:104
void build_it(SamplePointContainerParameters *sample_point_container_parameters_pt)
Helper function to actually build the thing.
Definition: mesh_as_geometric_object.h:96
void create_sample_point_container_parameters(Mesh *mesh_pt, SamplePointContainerParameters *&sample_point_container_parameters_pt)
Definition: mesh_as_geometric_object.cc:64

References build_it(), and oomph::MeshAsGeomObject_Helper::create_sample_point_container_parameters().

◆ MeshAsGeomObject() [2/4]

oomph::MeshAsGeomObject::MeshAsGeomObject ( SamplePointContainerParameters sample_point_container_parameters_pt)
inline

Constructor.

330  : GeomObject()
331  {
332  build_it(sample_point_container_parameters_pt);
333  }

References build_it().

◆ MeshAsGeomObject() [3/4]

oomph::MeshAsGeomObject::MeshAsGeomObject ( )
inline

Empty Constructor.

336 {}

◆ ~MeshAsGeomObject()

oomph::MeshAsGeomObject::~MeshAsGeomObject ( )
inline

Destructor.

340  {
342  }
SamplePointContainer * Sample_point_container_pt
Pointer to the sample point container.
Definition: mesh_as_geometric_object.h:270

References Sample_point_container_pt.

◆ MeshAsGeomObject() [4/4]

oomph::MeshAsGeomObject::MeshAsGeomObject ( const MeshAsGeomObject )
delete

Broken copy constructor.

Member Function Documentation

◆ build_it()

void oomph::MeshAsGeomObject::build_it ( SamplePointContainerParameters sample_point_container_parameters_pt)
inlineprivate

Helper function to actually build the thing.

98  {
99  Mesh_pt = sample_point_container_parameters_pt->mesh_pt();
100  if (dynamic_cast<RefineableBinArrayParameters*>(
101  sample_point_container_parameters_pt) != 0)
102  {
104  }
105  else if (dynamic_cast<NonRefineableBinArrayParameters*>(
106  sample_point_container_parameters_pt) != 0)
107  {
109  }
110 #ifdef OOMPH_HAS_CGAL
111  else if (dynamic_cast<CGALSamplePointContainerParameters*>(
112  sample_point_container_parameters_pt) != 0)
113  {
114  Sample_point_container_version = UseCGALSamplePointContainer;
115  }
116 #endif
117  else
118  {
119  throw OomphLibError("Wrong sample_point_container_parameters_pt",
122  }
123 
124 #ifdef OOMPH_HAS_MPI
125 
126  // Set communicator
127  Communicator_pt = Mesh_pt->communicator_pt();
128 
129 #endif
130 
131 
132  // Storage for the Lagrangian and Eulerian dimension
133  int dim[2] = {0, 0};
134 
135  // Set the Lagrangian dimension from the dimension of the first element
136  // if it exists (if not the Lagrangian dimension will be zero)
137  if (Mesh_pt->nelement() != 0)
138  {
139  dim[0] = Mesh_pt->finite_element_pt(0)->dim();
140  }
141 
142  // Read out the Eulerian dimension from the first node, if it exists.
143  //(if not the Eulerian dimension will be zero);
144  if (Mesh_pt->nnode() != 0)
145  {
146  dim[1] = Mesh_pt->node_pt(0)->ndim();
147  }
148 
149  // Need to do an Allreduce to ensure that the dimension is consistent
150  // even when no elements are assigned to a certain processor
151 #ifdef OOMPH_HAS_MPI
152 
153  // Only a problem if the mesh has been distributed
155  {
156  // Need a non-null communicator
157  if (Communicator_pt != 0)
158  {
159  int n_proc = Communicator_pt->nproc();
160  if (n_proc > 1)
161  {
162  int dim_reduce[2];
163  MPI_Allreduce(&dim,
164  &dim_reduce,
165  2,
166  MPI_INT,
167  MPI_MAX,
168  Communicator_pt->mpi_comm());
169 
170  dim[0] = dim_reduce[0];
171  dim[1] = dim_reduce[1];
172  }
173  }
174  }
175 #endif
176 
177  // Set the Lagrangian and Eulerian dimensions within this geometric object
178  this->set_nlagrangian_and_ndim(static_cast<unsigned>(dim[0]),
179  static_cast<unsigned>(dim[1]));
180 
181  // Create temporary storage for geometric Data (don't count
182  // Data twice!
183  std::set<Data*> tmp_geom_data;
184 
185  // Copy all the elements in the mesh into local storage
186  // N.B. elements must be able to have a geometric object representation.
187  unsigned n_sub_object = Mesh_pt->nelement();
188  Sub_geom_object_pt.resize(n_sub_object);
189  for (unsigned e = 0; e < n_sub_object; e++)
190  {
191  // (Try to) cast to a finite element:
193  dynamic_cast<FiniteElement*>(Mesh_pt->element_pt(e));
194 
195 #ifdef PARANOID
196  if (Sub_geom_object_pt[e] == 0)
197  {
198  std::ostringstream error_message;
199  error_message << "Unable to dynamic cast element: " << std::endl
200  << "into a FiniteElement: GeomObject representation is "
201  "not possible\n";
202  throw OomphLibError(error_message.str(),
205  }
206 #endif
207 
208  // Add the geometric Data of each element into set
209  unsigned ngeom_data = Sub_geom_object_pt[e]->ngeom_data();
210  for (unsigned i = 0; i < ngeom_data; i++)
211  {
212  tmp_geom_data.insert(Sub_geom_object_pt[e]->geom_data_pt(i));
213  }
214  }
215 
216  // Now copy unique geom Data values across into vector
217  unsigned ngeom = tmp_geom_data.size();
218  Geom_data_pt.resize(ngeom);
219  typedef std::set<Data*>::iterator IT;
220  unsigned count = 0;
221  for (IT it = tmp_geom_data.begin(); it != tmp_geom_data.end(); it++)
222  {
223  Geom_data_pt[count] = *it;
224  count++;
225  }
226 
227  // Build the right type of bin array
229  {
231 
233  new RefineableBinArray(sample_point_container_parameters_pt);
234  break;
235 
237 
239  new NonRefineableBinArray(sample_point_container_parameters_pt);
240  break;
241 
242 #ifdef OOMPH_HAS_CGAL
243 
244  case UseCGALSamplePointContainer:
245 
247  new CGALSamplePointContainer(sample_point_container_parameters_pt);
248  break;
249 
250 #endif
251 
252  default:
253 
254  oomph_info << "Sample_point_container_version = "
255  << Sample_point_container_version << std::endl;
256  throw OomphLibError("Sample_point_container_version",
259  }
260  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
NonRefineableBinArray class.
Definition: sample_point_container.h:819
RefineableBinArray class.
Definition: sample_point_container.h:521
unsigned dim() const
Definition: elements.h:2611
void set_nlagrangian_and_ndim(const unsigned &n_lagrangian, const unsigned &n_dim)
Set # of Lagrangian and Eulerian coordinates.
Definition: geom_objects.h:183
unsigned Sample_point_container_version
Definition: mesh_as_geometric_object.h:284
unsigned ngeom_data() const
How many items of Data does the shape of the object depend on?
Definition: mesh_as_geometric_object.h:351
Vector< FiniteElement * > Sub_geom_object_pt
Internal storage for the elements that constitute the object.
Definition: mesh_as_geometric_object.h:267
Data * geom_data_pt(const unsigned &j)
Definition: mesh_as_geometric_object.h:358
Mesh * Mesh_pt
Pointer to mesh.
Definition: mesh_as_geometric_object.h:280
Vector< Data * > Geom_data_pt
Vector of pointers to Data items that affects the object's shape.
Definition: mesh_as_geometric_object.h:264
bool is_mesh_distributed() const
Boolean to indicate if Mesh has been distributed.
Definition: mesh.h:1588
OomphCommunicator * communicator_pt() const
Definition: mesh.h:1600
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
Node *& node_pt(const unsigned long &n)
Return pointer to global node n.
Definition: mesh.h:436
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
@ UseRefineableBinArray
Definition: sample_point_parameters.h:41
@ UseNonRefineableBinArray
Definition: sample_point_parameters.h:42
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::Mesh::communicator_pt(), oomph::FiniteElement::dim(), e(), oomph::Mesh::element_pt(), oomph::Mesh::finite_element_pt(), Geom_data_pt, geom_data_pt(), i, oomph::Mesh::is_mesh_distributed(), Mesh_pt, oomph::SamplePointContainerParameters::mesh_pt(), oomph::Node::ndim(), oomph::Mesh::nelement(), ngeom_data(), oomph::Mesh::nnode(), oomph::Mesh::node_pt(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, Sample_point_container_pt, Sample_point_container_version, oomph::GeomObject::set_nlagrangian_and_ndim(), Sub_geom_object_pt, oomph::UseNonRefineableBinArray, and oomph::UseRefineableBinArray.

Referenced by MeshAsGeomObject().

◆ dposition()

void oomph::MeshAsGeomObject::dposition ( const Vector< double > &  xi,
DenseMatrix< double > &  drdxi 
) const
inlinevirtual

Return the derivative of the position.

Reimplemented from oomph::GeomObject.

444  {
445  throw OomphLibError("dposition() not implemented",
448  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ finite_element_pt()

FiniteElement* oomph::MeshAsGeomObject::finite_element_pt ( const unsigned e)
inline

Return pointer to e-th finite element.

295  {
296  return Sub_geom_object_pt[e];
297  }

References e(), and Sub_geom_object_pt.

◆ geom_data_pt()

Data* oomph::MeshAsGeomObject::geom_data_pt ( const unsigned j)
inlinevirtual

Return pointer to the j-th Data item that the object's shape depends on

Reimplemented from oomph::GeomObject.

359  {
360  return Geom_data_pt[j];
361  }
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Geom_data_pt, and j.

Referenced by build_it().

◆ locate_zeta()

void oomph::MeshAsGeomObject::locate_zeta ( const Vector< double > &  zeta,
GeomObject *&  sub_geom_object_pt,
Vector< double > &  s,
const bool use_coordinate_as_initial_guess = false 
)
inlinevirtual

Find the sub geometric object and local coordinate therein that corresponds to the intrinsic coordinate zeta. If sub_geom_object_pt=0 on return from this function, none of the constituent sub-objects contain the required coordinate. Following from the general interface to this function in GeomObjects, setting the optional bool argument to true means that each time the sub-object's locate_zeta function is called, the coordinate argument "s" is used as the initial guess. However, this doesn't make sense here and the argument is ignored (though a warning is issued when the code is compiled in PARANOID setting)

Reimplemented from oomph::GeomObject.

377  {
378 #ifdef PARANOID
379  if (use_coordinate_as_initial_guess)
380  {
381  OomphLibWarning(
382  "Ignoring the use_coordinate_as_initial_guess argument.",
383  "MeshAsGeomObject::locate_zeta()",
385  }
386 #endif
387 
388 
389  // Do locate in sample point container
390  Sample_point_container_pt->locate_zeta(zeta, sub_geom_object_pt, s);
391  }
virtual void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s)=0
RealScalar s
Definition: level1_cplx_impl.h:130
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 SamplePointContainer::locate_zeta(), OOMPH_EXCEPTION_LOCATION, s, Sample_point_container_pt, and Eigen::zeta().

Referenced by DarcyProblem< ELEMENT >::doc_shape_functions(), position(), SegregatedFSICollapsibleChannelProblem< ELEMENT >::SegregatedFSICollapsibleChannelProblem(), SegregatedFSIDrivenCavityProblem< ELEMENT >::SegregatedFSIDrivenCavityProblem(), oomph::LineVisualiser::setup(), oomph::HelmholtzMGPreconditioner< DIM >::setup_interpolation_matrices(), oomph::MGSolver< DIM >::setup_interpolation_matrices_unstructured(), and oomph::HelmholtzMGPreconditioner< DIM >::setup_interpolation_matrices_unstructured().

◆ nelement()

unsigned oomph::MeshAsGeomObject::nelement ( )
inline

Number of elements in the underlying mesh.

309  {
310  return Sub_geom_object_pt.size();
311  }

References Sub_geom_object_pt.

◆ ngeom_data()

unsigned oomph::MeshAsGeomObject::ngeom_data ( ) const
inlinevirtual

How many items of Data does the shape of the object depend on?

Reimplemented from oomph::GeomObject.

352  {
353  return Geom_data_pt.size();
354  }

References Geom_data_pt.

Referenced by build_it().

◆ operator=()

void oomph::MeshAsGeomObject::operator= ( const MeshAsGeomObject )
delete

Broken assignment operator.

◆ position() [1/2]

void oomph::MeshAsGeomObject::position ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  r 
) const
inlinevirtual

Parametrised position on object: r(zeta). Evaluated at previous timestep. t=0: current time; t>0: previous timestep. This provides an (expensive!) default implementation in which we loop over all the constituent sub-objects and check if they contain zeta and then evaluate their position() function.

Reimplemented from oomph::GeomObject.

412  {
413  // Storage for the GeomObject that contains the zeta coordinate
414  // and the intrinsic coordinate within it.
415  GeomObject* sub_geom_object_pt;
416  const unsigned n_lagrangian = this->nlagrangian();
417  Vector<double> s(n_lagrangian);
418 
419  // Find the sub object containing zeta, and the local intrinsic coordinate
420  // within it
421  const_cast<MeshAsGeomObject*>(this)->locate_zeta(
422  zeta, sub_geom_object_pt, s);
423  if (sub_geom_object_pt == 0)
424  {
425  std::ostringstream error_message;
426  error_message << "Cannot locate zeta ";
427  for (unsigned i = 0; i < n_lagrangian; i++)
428  {
429  error_message << zeta[i] << " ";
430  }
431  error_message << std::endl;
432  Mesh_pt->output("most_recent_mesh.dat");
433  throw OomphLibError(error_message.str(),
436  }
437  // Call that sub-object's position function
438  sub_geom_object_pt->position(t, s, r);
439 
440  } // end of position
unsigned nlagrangian() const
Access function to # of Lagrangian coordinates.
Definition: geom_objects.h:171
MeshAsGeomObject()
Empty Constructor.
Definition: mesh_as_geometric_object.h:336
void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Definition: mesh_as_geometric_object.h:373
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
r
Definition: UniformPSDSelfTest.py:20
t
Definition: plotPSD.py:36

References i, locate_zeta(), Mesh_pt, oomph::GeomObject::nlagrangian(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Mesh::output(), oomph::GeomObject::position(), UniformPSDSelfTest::r, s, plotPSD::t, and Eigen::zeta().

◆ position() [2/2]

void oomph::MeshAsGeomObject::position ( const Vector< double > &  zeta,
Vector< double > &  r 
) const
inlinevirtual

Return the position as a function of the intrinsic coordinate zeta. This provides an (expensive!) default implementation in which we loop over all the constituent sub-objects and check if they contain zeta and then evaluate their position() function.

Implements oomph::GeomObject.

398  {
399  // Call position function at current timestep:
400  unsigned t = 0;
401  position(t, zeta, r);
402  }
void position(const Vector< double > &zeta, Vector< double > &r) const
Definition: mesh_as_geometric_object.h:397

References UniformPSDSelfTest::r, plotPSD::t, and Eigen::zeta().

Referenced by oomph::RefineableTetgenMesh< ELEMENT >::snap_nodes_onto_boundary().

◆ sample_point_container_pt()

SamplePointContainer* oomph::MeshAsGeomObject::sample_point_container_pt ( ) const
inline

Pointer to the sample point container.

289  {
291  }

References Sample_point_container_pt.

Referenced by oomph::RefineableGmshTetMesh< ELEMENT >::adapt(), and oomph::LineVisualiser::setup().

◆ sample_point_container_version()

unsigned oomph::MeshAsGeomObject::sample_point_container_version ( ) const
inline

Which sample point container is used in locate zeta? (uses enum Sample_Point_Container_Type)

303  {
305  }

References Sample_point_container_version.

Member Data Documentation

◆ Geom_data_pt

Vector<Data*> oomph::MeshAsGeomObject::Geom_data_pt
private

Vector of pointers to Data items that affects the object's shape.

Referenced by build_it(), geom_data_pt(), and ngeom_data().

◆ Mesh_pt

Mesh* oomph::MeshAsGeomObject::Mesh_pt
private

Pointer to mesh.

Referenced by build_it(), and position().

◆ Sample_point_container_pt

SamplePointContainer* oomph::MeshAsGeomObject::Sample_point_container_pt
private

Pointer to the sample point container.

Referenced by build_it(), locate_zeta(), sample_point_container_pt(), and ~MeshAsGeomObject().

◆ Sample_point_container_version

unsigned oomph::MeshAsGeomObject::Sample_point_container_version
private

Which version of the sample point container are we using?

Referenced by build_it(), and sample_point_container_version().

◆ Sub_geom_object_pt

Vector<FiniteElement*> oomph::MeshAsGeomObject::Sub_geom_object_pt
private

Internal storage for the elements that constitute the object.

Referenced by build_it(), finite_element_pt(), and nelement().


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