oomph::GeomObject Class Referenceabstract

#include <geom_objects.h>

+ Inheritance diagram for oomph::GeomObject:

Public Member Functions

 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 unsigned ngeom_data () const
 
virtual Datageom_data_pt (const unsigned &j)
 
virtual void position (const Vector< double > &zeta, Vector< double > &r) const =0
 Parametrised position on object at current time: r(zeta). More...
 
virtual void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) 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 dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
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 locate_zeta (const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 

Protected Attributes

unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Detailed Description

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// A geometric object is an object that provides a parametrised description of its shape via the function GeomObject::position(...).

The minimum functionality is: The geometric object has a number of Lagrangian (intrinsic) coordinates that parametrise the (Eulerian) position vector, whose dimension might differ from the number of Lagrangian (intrinsic) coordinates (e.g. for shell-like objects).

We might also need the derivatives of the position Vector w.r.t. to the Lagrangian (intrinsic) coordinates and interfaces for this functionality are provided. [Note: For some geometric objects it might be too tedious to work out the derivatives and they might not be needed anyway. In other cases we might always need the position vector and all derivatives at the same time. We provide suitable interfaces for these cases in virtual but broken (rather than pure virtual) form so the user can (but doesn't have to) provide the required versions by overloading.]

The shape of a geometric object is usually determined by a number of parameters whose value might have to be determined as part of the overall solution (e.g. for geometric objects that represent elastic walls). The geometric object therefore has a vector of (pointers to) geometric Data, which can be free/pinned and have a time history, etc. This makes it possible to ‘upgrade’ GeomObjects to GeneralisedElements – in this case the geometric Data plays the role of internal Data in the GeneralisedElement. Conversely, FiniteElements, in which a geometry (spatial coordinate) has been defined, inherit from GeomObjects, which is particularly useful in FSI computations: Meshing of moving domains is typically performed by representing the domain as an object of type Domain and, by default, Domain boundaries are represented by GeomObjects. In FSI computations, the boundary of the fluid domain is represented by a number of solid mechanics elements. These elements are, in fact, GeomObjects via inheritance so that the we can use the standard interfaces of the GeomObject class for mesh generation. An example is the class FSIHermiteBeamElement which is derived from the class HermiteBeamElement (a ‘normal’ beam element) and the GeomObject class.

The shape of a geometric object can have an explicit time-dependence, for instance in cases where a domain boundary is performing prescribed motions. We provide access to the ‘global’ time by giving the object a pointer to a timestepping scheme. [Note that, within the overall FE code, time is only ever evaluated at discrete instants (which are accessible via the timestepper), never in continuous form]. The timestepper is also needed to evaluate time-derivatives if the geometric Data carries a time history.

Constructor & Destructor Documentation

◆ GeomObject() [1/5]

oomph::GeomObject::GeomObject ( )
inline

Default constructor.

unsigned NLagrangian
Number of Lagrangian (intrinsic) coordinates.
Definition: geom_objects.h:428
TimeStepper * Geom_object_time_stepper_pt
Definition: geom_objects.h:435
unsigned Ndim
Number of Eulerian coordinates.
Definition: geom_objects.h:431

◆ GeomObject() [2/5]

oomph::GeomObject::GeomObject ( const unsigned ndim)
inline

Constructor: Pass dimension of geometric object (# of Eulerian coords = # of Lagrangian coords; no time history available/needed)

110  {
111  }
unsigned ndim() const
Access function to # of Eulerian coordinates.
Definition: geom_objects.h:177

◆ GeomObject() [3/5]

oomph::GeomObject::GeomObject ( const unsigned nlagrangian,
const unsigned ndim 
)
inline

Constructor: pass # of Eulerian and Lagrangian coordinates. No time history available/needed

118  {
119 #ifdef PARANOID
120  if (nlagrangian > ndim)
121  {
122  std::ostringstream error_message;
123  error_message << "# of Lagrangian coordinates " << nlagrangian
124  << " cannot be bigger than # of Eulerian ones " << ndim
125  << std::endl;
126 
127  throw OomphLibError(error_message.str(),
130  }
131 #endif
132  }
unsigned nlagrangian() const
Access function to # of Lagrangian coordinates.
Definition: geom_objects.h:171
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References ndim(), nlagrangian(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ GeomObject() [4/5]

oomph::GeomObject::GeomObject ( const unsigned nlagrangian,
const unsigned ndim,
TimeStepper time_stepper_pt 
)
inline

Constructor: pass # of Eulerian and Lagrangian coordinates and pointer to time-stepper which is used to handle the position at previous timesteps and allows the evaluation of veloc/acceleration etc. in cases where the GeomData varies with time.

143  Ndim(ndim),
145  {
146 #ifdef PARANOID
147  if (nlagrangian > ndim)
148  {
149  std::ostringstream error_message;
150  error_message << "# of Lagrangian coordinates " << nlagrangian
151  << " cannot be bigger than # of Eulerian ones " << ndim
152  << std::endl;
153 
154  throw OomphLibError(error_message.str(),
157  }
158 #endif
159  }
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192

References ndim(), nlagrangian(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ GeomObject() [5/5]

oomph::GeomObject::GeomObject ( const GeomObject dummy)
delete

Broken copy constructor.

◆ ~GeomObject()

virtual oomph::GeomObject::~GeomObject ( )
inlinevirtual

(Empty) destructor

168 {}

Member Function Documentation

◆ d2position() [1/2]

virtual void oomph::GeomObject::d2position ( const Vector< double > &  zeta,
RankThreeTensor< double > &  ddrdzeta 
) const
inlinevirtual

◆ d2position() [2/2]

virtual void oomph::GeomObject::d2position ( const Vector< double > &  zeta,
Vector< double > &  r,
DenseMatrix< double > &  drdzeta,
RankThreeTensor< double > &  ddrdzeta 
) const
inlinevirtual

Posn Vector and its 1st & 2nd derivatives w.r.t. to coordinates: \( \frac{dR_i}{d \zeta_\alpha}\) = drdzeta(alpha,i). \( \frac{d^2R_i}{d \zeta_\alpha d \zeta_\beta}\) = ddrdzeta(alpha,beta,i). Evaluated at current time.

Reimplemented in oomph::EllipticalTube, oomph::Ellipse, oomph::StraightLine, FlatPlate, UndeformedWall, SpikedLine, SinusoidalWall, UndeformedLeaflet, UndeformedWall, UndeformedLeaflet, UndeformedLeaflet, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedWall, and UndeformedLeaflet.

361  {
362  throw OomphLibError(
363  "You must specify d2position() for your own object! \n",
366  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ dposition()

virtual void oomph::GeomObject::dposition ( const Vector< double > &  zeta,
DenseMatrix< double > &  drdzeta 
) const
inlinevirtual

Derivative of position Vector w.r.t. to coordinates: \( \frac{dR_i}{d \zeta_\alpha}\) = drdzeta(alpha,i). Evaluated at current time.

Reimplemented in oomph::Ellipse, oomph::StraightLine, SpikedLine, SinusoidalWall, and oomph::MeshAsGeomObject.

329  {
330  throw OomphLibError(
331  "You must specify dposition() for your own object! \n",
334  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ dposition_dt()

virtual void oomph::GeomObject::dposition_dt ( const Vector< double > &  zeta,
const unsigned j,
Vector< double > &  drdt 
)
inlinevirtual

j-th time-derivative on object at current time: \( \frac{d^{j} r(\zeta)}{dt^j} \).

Reimplemented in oomph::FiniteElement, oomph::ImmersedRigidBodyElement, oomph::PseudoBucklingRing, oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >, and AxisymOscillatingDisk.

295  {
296  // If the index is zero the return the position
297  if (j == 0)
298  {
299  position(zeta, drdt);
300  }
301  // Otherwise assume that the geometric object is static
302  // and return zero after throwing a warning
303  else
304  {
305  std::ostringstream warning_stream;
306  warning_stream
307  << "Using default (static) assignment " << j
308  << "-th time derivative in GeomObject::dposition_dt(...) is zero\n"
309  << "Overload for your specific geometric object if this is not \n"
310  << "appropriate. \n";
311  OomphLibWarning(warning_stream.str(),
312  "GeomObject::dposition_dt()",
314 
315  unsigned n = drdt.size();
316  for (unsigned i = 0; i < n; i++)
317  {
318  drdt[i] = 0.0;
319  }
320  }
321  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
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
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, j, n, OOMPH_EXCEPTION_LOCATION, position(), and Eigen::zeta().

Referenced by oomph::SolidFiniteElement::fill_in_generic_jacobian_for_solid_ic().

◆ geom_data_pt()

virtual Data* oomph::GeomObject::geom_data_pt ( const unsigned j)
inlinevirtual

Return pointer to the j-th Data item that the object's shape depends on. This is implemented as a broken virtual function. You must overload this for GeomObjects that contain geometric Data, i.e. GeomObjects whose shape depends on Data that may contain unknowns in the overall Problem.

Reimplemented in oomph::ImmersedRigidBodyElement, oomph::RefineableSolidElement, oomph::PseudoBucklingRing, oomph::MeshAsGeomObject, oomph::Circle, oomph::Ellipse, oomph::StraightLine, oomph::FaceElementAsGeomObject< ELEMENT >, oomph::SolidFiniteElement, oomph::FiniteElement, oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >, SpikedLine, SinusoidalWall, oomph::GeneralCircle, and GeneralCircle.

234  {
235  std::ostringstream error_message;
236  error_message
237  << "GeomObject::geom_data_pt() is a broken virtual function.\n"
238  << "Please implement it (and its companion GeomObject::ngeom_data())\n"
239  << "for any GeomObject whose shape depends on Data whose values may \n"
240  << "be unknowns in the global Problem. \n"
241  << "If you have arrived here in a parallel job then it may be the case "
242  "\n"
243  << "that you have not set the keep_all_elements_as_halos() flag to "
244  "true \n"
245  << "for the MeshAsGeomObject representing the lower-dimensional mesh \n"
246  << "in a problem with multiple meshes. \n";
247  throw OomphLibError(
248  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
249  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data(), AirwayReopeningProblem< ELEMENT >::connect_walls(), and oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt().

◆ interpolated_zeta()

virtual void oomph::GeomObject::interpolated_zeta ( const Vector< double > &  s,
Vector< double > &  zeta 
) const
inlinevirtual

A geometric object may be composed of many sub-objects each with their own local coordinate. This function returns the "global" intrinsic coordinate zeta (within the compound object), at a given local coordinate s (i.e. the intrinsic coordinate of the sub-GeomObject. In simple (non-compound) GeomObjects, the local intrinsic coordinate is the global intrinsic coordinate and so the function merely returns s. To make it less likely that the default implementation is called in error (because it is not overloaded in a derived GeomObject where the default is not appropriate, we do at least check that s and zeta have the same size if called in PARANOID mode.

Reimplemented in oomph::FiniteElement.

406  {
407 #ifdef PARANOID
408  if (zeta.size() != s.size())
409  {
410  std::ostringstream error_message;
411  error_message << "You've called the default implementation of "
412  << "GeomObject::interpolated_zeta() \n"
413  << "but zeta.size()=" << zeta.size()
414  << "and s.size()=" << s.size() << std::endl
415  << "This doesn't make sense! You probably have to \n"
416  << "overload this function in your specific GeomObject\n";
417  throw OomphLibError(error_message.str(),
420  }
421 #endif
422  // By default the global intrinsic coordinate is equal to the local one
423  zeta = s;
424  }
RealScalar s
Definition: level1_cplx_impl.h:130

References OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, s, and Eigen::zeta().

◆ locate_zeta()

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

A geometric object may be composed of may sub-objects (e.g. a finite-element representation of a boundary). In order to implement sparse update functions, it is necessary to know the sub-object and local coordinate within that sub-object at a given intrinsic coordinate, zeta. Note that only one sub-object can "cover" any given intrinsic position. If the position is at an "interface" between sub-objects, either one can be returned. The default implementation merely returns, the pointer to the "entire" GeomObject and the coordinate, zeta The optional boolean flag only applies if a Newton method is used to find the value of zeta, and if true the value of the coordinate s is used as the initial guess for the method. If the flag is false (the default) a value of s=0 is used as the initial guess.

Reimplemented in oomph::MeshAsGeomObject, oomph::FiniteElement, and oomph::FSIHermiteBeamElement.

386  {
387  // By default, the local coordinate is intrinsic coordinate
388  s = zeta;
389  // The sub_object is the entire object
390  sub_geom_object_pt = this;
391  }

References s, and Eigen::zeta().

Referenced by SCoupling< M, O >::computeSCouplingForcesFromTriangles(), AirwayReopeningProblem< ELEMENT >::connect_walls(), SolidProblem< ELEMENT_TYPE >::get_x(), PressureWaveFSIProblem< FLUID_ELEMENT, SOLID_ELEMENT >::PressureWaveFSIProblem(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::setup_algebraic_node_update(), oomph::RefineableAlgebraicChannelWithLeafletMesh< ELEMENT >::update_node_update(), oomph::RefineableAlgebraicCollapsibleChannelMesh< ELEMENT >::update_node_update(), oomph::RefineableAlgebraicCylinderWithFlagMesh< ELEMENT >::update_node_update(), oomph::AlgebraicFishMesh< ELEMENT >::update_node_update(), and oomph::RefineableAlgebraicFSIDrivenCavityMesh< ELEMENT >::update_node_update().

◆ ndim()

unsigned oomph::GeomObject::ndim ( ) const
inline

Access function to # of Eulerian coordinates.

178  {
179  return Ndim;
180  }

References Ndim.

Referenced by oomph::FSIHermiteBeamElement::dposition_dlagrangian_at_local_coordinate(), oomph::FSIDiagHermiteShellElement::dposition_dlagrangian_at_local_coordinate(), oomph::NavierStokesSurfacePowerElement< ELEMENT >::drag_force(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), GeomObject(), oomph::NavierStokesSurfaceDragTorqueElement< ELEMENT >::get_drag_and_torque(), oomph::KirchhoffLoveShellEquations::get_energy(), oomph::KirchhoffLoveBeamEquations::get_energy(), oomph::NavierStokesSurfacePowerElement< ELEMENT >::get_kinetic_energy_flux(), oomph::KirchhoffLoveBeamEquations::get_non_unit_tangent(), oomph::KirchhoffLoveBeamEquations::get_normal(), oomph::NavierStokesSurfacePowerElement< ELEMENT >::get_rate_of_traction_work(), oomph::NavierStokesSurfacePowerElement< ELEMENT >::get_rate_of_traction_work_components(), oomph::NavierStokesSurfacePowerElement< ELEMENT >::get_volume_flux(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_central_region(), oomph::AlgebraicFishMesh< ELEMENT >::node_update_in_body(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_central_box(), oomph::AlgebraicFishMesh< ELEMENT >::node_update_in_fin(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_lower_right_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_upper_left_box(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_lower_right_region(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_upper_left_region(), oomph::HermiteBeamElement::output(), oomph::PMLTimeHarmonicLinearElasticityTractionElement< ELEMENT >::output(), oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::output(), and oomph::DiskLikeGeomObjectWithBoundaries::zeta_on_boundary().

◆ ngeom_data()

virtual unsigned oomph::GeomObject::ngeom_data ( ) const
inlinevirtual

How many items of Data does the shape of the object depend on? This is implemented as a broken virtual function. You must overload this for GeomObjects that contain geometric Data, i.e. GeomObjects whose shape depends on Data that may contain unknowns in the overall Problem.

Reimplemented in oomph::ImmersedRigidBodyElement, oomph::RefineableSolidElement, oomph::PseudoBucklingRing, oomph::MeshAsGeomObject, oomph::EllipticalTube, oomph::Circle, oomph::Ellipse, oomph::StraightLine, oomph::FaceElementAsGeomObject< ELEMENT >, oomph::SolidFiniteElement, oomph::FiniteElement, oomph::ElementWithMovingNodes, oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >, OscillatingWall, WarpedPlane, OscillatingWall, WarpedLine, WarpedLine, FlatPlate, Flag_definition::TipOfFlag, Flag_definition::BottomOfFlag, Flag_definition::TopOfFlag, SpikedLine, SinusoidalWall, OscillatingWall, OscillatingWall, OscillatingWall, Leaflet, UndeformedLeaflet, UndeformedLeaflet, WarpedPlane, WarpedLine, WarpedLine, UndeformedLeaflet, UndeformedLeaflet, oomph::GeneralCircle, WarpedLine, WarpedLine, and GeneralCircle.

210  {
211  std::ostringstream error_message;
212  error_message
213  << "GeomObject::ngeom_data() is a broken virtual function.\n"
214  << "Please implement it (and its companion "
215  "GeomObject::geom_data_pt())\n"
216  << "for any GeomObject whose shape depends on Data whose values may \n"
217  << "be unknowns in the global Problem. \n"
218  << "If you have arrived here in a parallel job then it may be the case "
219  "\n"
220  << "that you have not set the keep_all_elements_as_halos() flag to "
221  "true \n"
222  << "for the MeshAsGeomObject representing the lower-dimensional mesh \n"
223  << "in a problem with multiple meshes. \n";
224  throw OomphLibError(
225  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
226  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data(), AirwayReopeningProblem< ELEMENT >::connect_walls(), and oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt().

◆ nlagrangian()

◆ operator=()

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

Broken assignment operator.

◆ position() [1/3]

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

Parametrised position on object: r(zeta). Evaluated at the continuous time value, t.

Reimplemented in OscillatingCylinder, and OscillatingCylinder.

279  {
280  std::ostringstream error_message;
281  error_message << "GeomObject::position() is a broken virtual function.\n"
282  << "Please implement it for any GeomObject whose shape\n"
283  << "is time-dependent and will be used in the extrusion\n"
284  << "of a mesh (in the time direction).\n";
285  throw OomphLibError(
286  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
287  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ position() [2/3]

virtual void oomph::GeomObject::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. Works for t=0 but needs to be overloaded if genuine time-dependence is required.

Reimplemented in oomph::PseudoBucklingRing, oomph::MeshAsGeomObject, oomph::EllipticalTube, oomph::Circle, oomph::Ellipse, oomph::StraightLine, oomph::WarpedCircularDisk, oomph::FiniteElement, oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >, OscillatingWall, WarpedPlane, OscillatingWall, WarpedLine, WarpedLine, WarpedLine, FlatPlate, UndeformedWall, SpikedLine, SinusoidalWall, OscillatingWall, OscillatingWall, OscillatingWall, UndeformedLeaflet, UndeformedWall, UndeformedLeaflet, WarpedPlane, WarpedLine, WarpedLine, UndeformedLeaflet, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedLeaflet, UndeformedWall, oomph::GeneralCircle, oomph::SimpleCircle, WarpedLine, WarpedLine, oomph::ImmersedRigidBodyTriangleMeshPolygon, oomph::ImmersedRigidBodyElement, MyEllipse, MyUnitCircle, MyEllipse, UnitCircle, MyUnitCircle, MyEllipse, Flag_definition::TipOfFlag, Flag_definition::BottomOfFlag, Flag_definition::TopOfFlag, MyCylinder, OscillatingCylinder, MyEllipse, GeneralEllipse, OscillatingCylinder, MyHelicalCylinder, MyCurvedCylinder, Leaflet, MyEllipse, OscillatingCylinder, GeneralCircle, GeneralCircle, and FilledCircle.

262  {
263  if (t != 0)
264  {
265  throw OomphLibError(
266  "Calling steady position() from discrete unsteady position()",
269  }
270  position(zeta, r);
271  }
r
Definition: UniformPSDSelfTest.py:20
t
Definition: plotPSD.py:36

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

◆ position() [3/3]

virtual void oomph::GeomObject::position ( const Vector< double > &  zeta,
Vector< double > &  r 
) const
pure virtual

Parametrised position on object at current time: r(zeta).

Implemented in oomph::PseudoBucklingRing, oomph::MeshAsGeomObject, oomph::EllipticalTube, oomph::Circle, oomph::Ellipse, oomph::StraightLine, oomph::WarpedCircularDisk, oomph::FiniteElement, oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >, OscillatingWall, MyStraightLine, WarpedPlane, OscillatingWall, WarpedLine, WarpedLine, WarpedLine, FlatPlate, UndeformedWall, SpikedLine, SinusoidalWall, WavyWall, OscillatingWall, OscillatingWall, OscillatingWall, UndeformedLeaflet, UndeformedWall, UndeformedLeaflet, WarpedPlane, WarpedLine, WarpedLine, UndeformedLeaflet, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedWall, UndeformedLeaflet, UndeformedWall, oomph::GeneralCircle, oomph::SimpleCircle, MyStraightLine, MyStraightLine, WarpedLine, WarpedLine, oomph::ImmersedRigidBodyTriangleMeshPolygon, oomph::ImmersedRigidBodyElement, MyEllipse, MyUnitCircle, MyEllipse, AxisymOscillatingDisk, UnitCircle, MyUnitCircle, MyEllipse, Flag_definition::TipOfFlag, Flag_definition::BottomOfFlag, Flag_definition::TopOfFlag, MyCylinder, OscillatingCylinder, MyEllipse, GeneralEllipse, OscillatingCylinder, GeneralEllipse, MyHelicalCylinder, MyCurvedCylinder, Leaflet, GeneralEllipse, oomph::HalfEllipse, oomph::HalfEllipse, MyEllipse, OscillatingCylinder, GeneralEllipse, GeneralEllipse, GeneralCircle, GeneralCircle, FilledCircle, and oomph::YoungLaplaceEquations.

Referenced by oomph::MyAlgebraicCollapsibleChannelMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::AlgebraicChannelWithLeafletMesh(), oomph::ChannelWithLeafletDomain::ChannelWithLeafletDomain(), oomph::TriangleMeshCurveSection::connect_final_vertex_to_curviline(), oomph::TriangleMeshCurveSection::connect_initial_vertex_to_curviline(), oomph::DiskTetMeshFacetedSurface::DiskTetMeshFacetedSurface(), UnstructuredFluidProblem< ELEMENT >::doc_boundary_coordinates(), dposition_dt(), oomph::ImmersedRigidBodyElement::dposition_dt(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::TriangleMeshCurviLine::final_vertex_coordinate(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::find_distance_to_free_surface(), oomph::BrethertonSpineMesh< ELEMENT, SpineLineFluidInterfaceElement< ELEMENT > >::initial_element_reorder(), oomph::TriangleMeshCurviLine::initial_vertex_coordinate(), oomph::ChannelWithLeafletDomain::macro_bound_I_E(), oomph::ChannelWithLeafletDomain::macro_bound_I_W(), oomph::ChannelWithLeafletDomain::macro_bound_II_E(), oomph::ChannelWithLeafletDomain::macro_bound_II_W(), oomph::FullCircleDomain::macro_element_boundary(), oomph::TubeDomain::macro_element_boundary(), oomph::CylinderWithFlagDomain::macro_element_boundary(), oomph::RectangleWithHoleDomain::macro_element_boundary(), oomph::RectangleWithHoleAndAnnularRegionDomain::macro_element_boundary(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_central_region(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::node_update_I(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_I(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::node_update_II(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_II(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_III(), oomph::AlgebraicFishMesh< ELEMENT >::node_update_in_body(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_central_box(), oomph::AlgebraicFishMesh< ELEMENT >::node_update_in_fin(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_lower_right_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_upper_left_box(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_IV(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_IX(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_lower_right_region(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_upper_left_region(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_V(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_VI(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_VII(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_VIII(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::output(), oomph::TriangleMeshCurviLine::output(), position(), oomph::MeshAsGeomObject::position(), oomph::ImmersedRigidBodyElement::position(), oomph::DiskLikeGeomObjectWithBoundaries::position_on_boundary(), oomph::RectangleWithHoleAndAnnularRegionDomain::project_point_on_cylinder_to_annular_boundary(), oomph::QuarterPipeDomain::r_B(), oomph::QuarterTubeDomain::r_bot_right_B(), oomph::QuarterTubeDomain::r_bot_right_D(), oomph::QuarterCircleSectorDomain::r_bot_right_E(), oomph::QuarterTubeDomain::r_bot_right_F(), oomph::QuarterTubeDomain::r_bot_right_R(), oomph::QuarterCircleSectorDomain::r_bot_right_S(), oomph::QuarterTubeDomain::r_bot_right_U(), oomph::QuarterCircleSectorDomain::r_bot_right_W(), oomph::QuarterTubeDomain::r_centr_B(), oomph::QuarterTubeDomain::r_centr_D(), oomph::QuarterTubeDomain::r_centr_F(), oomph::QuarterTubeDomain::r_centr_L(), oomph::QuarterTubeDomain::r_centr_R(), oomph::QuarterCircleSectorDomain::r_centr_S(), oomph::QuarterTubeDomain::r_centr_U(), oomph::QuarterCircleSectorDomain::r_centr_W(), oomph::QuarterPipeDomain::r_D(), oomph::CollapsibleChannelDomain::r_E_collapsible(), oomph::QuarterPipeDomain::r_F(), oomph::QuarterPipeDomain::r_L(), oomph::CollapsibleChannelDomain::r_N_collapsible(), oomph::QuarterPipeDomain::r_R(), oomph::CollapsibleChannelDomain::r_S_collapsible(), oomph::QuarterTubeDomain::r_top_left_B(), oomph::QuarterCircleSectorDomain::r_top_left_E(), oomph::QuarterTubeDomain::r_top_left_F(), oomph::QuarterTubeDomain::r_top_left_L(), oomph::QuarterCircleSectorDomain::r_top_left_N(), oomph::QuarterCircleSectorDomain::r_top_left_S(), oomph::QuarterTubeDomain::r_top_left_U(), oomph::QuarterCircleSectorDomain::r_top_left_W(), oomph::QuarterPipeDomain::r_U(), oomph::FishDomain::r_upper_body_E(), oomph::FishDomain::r_upper_body_N(), oomph::FishDomain::r_upper_body_S(), oomph::FishDomain::r_upper_body_W(), oomph::FishDomain::r_upper_fin_N(), oomph::FishDomain::r_upper_fin_S(), oomph::FishDomain::r_upper_fin_W(), oomph::CollapsibleChannelDomain::r_W_collapsible(), oomph::RefineableImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::RefineableQuadMeshWithMovingCylinder< MyRefineableQTaylorHoodElement >::RefineableQuadMeshWithMovingCylinder(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::reposition_spines(), oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::setup_algebraic_node_update(), oomph::UnstructuredTwoDMeshGeometryBase::setup_boundary_coordinates(), oomph::ChannelWithLeafletDomain::slanted_bound_up(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), oomph::UnstructuredTwoDMeshGeometryBase::snap_nodes_onto_geometric_objects(), oomph::TetMeshBase::snap_to_quadratic_surface(), oomph::ChannelSpineMesh< ELEMENT >::spine_node_update(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_channel(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_film_lower(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_film_upper(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_horizontal_transition_lower(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_horizontal_transition_upper(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_vertical_transition_lower(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_vertical_transition_upper(), and oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::square_of_l2_norm_of_error().

◆ set_nlagrangian_and_ndim()

void oomph::GeomObject::set_nlagrangian_and_ndim ( const unsigned n_lagrangian,
const unsigned n_dim 
)
inline

Set # of Lagrangian and Eulerian coordinates.

185  {
186  NLagrangian = n_lagrangian;
187  Ndim = n_dim;
188  }

References Ndim, and NLagrangian.

Referenced by oomph::MeshAsGeomObject::build_it(), oomph::FaceElementAsGeomObject< ELEMENT >::FaceElementAsGeomObject(), and oomph::FSIWallElement::setup_fsi_wall_element().

◆ time_stepper_pt() [1/2]

TimeStepper*& oomph::GeomObject::time_stepper_pt ( )
inline

Access function for pointer to time stepper: Null if object is not time-dependent

193  {
195  }

References Geom_object_time_stepper_pt.

Referenced by oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), oomph::SurfactantTransportInterfaceElement::add_additional_residual_contributions_interface(), oomph::Circle::Circle(), oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::construct_boundary_node(), oomph::FiniteElement::construct_boundary_node(), oomph::SolidFiniteElement::construct_boundary_node(), oomph::DGElement::construct_boundary_nodes_and_faces(), oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::construct_node(), oomph::FiniteElement::construct_node(), oomph::SolidFiniteElement::construct_node(), oomph::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::PoroelasticityEquations< DIM >::d2u_dt2(), oomph::AxisymmetricLinearElasticityEquationsBase::d2u_dt2_axisymmetric_linear_elasticity(), oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave(), oomph::LinearElasticityEquationsBase< DIM >::d2u_dt2_linear_elasticity(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react(), oomph::SurfactantTransportInterfaceElement::dcdt_surface(), oomph::PerturbedSpineLinearisedAxisymmetricFluidInterfaceElement< ELEMENT >::dH_dt(), oomph::LinearisedAxisymmetricNavierStokesEquations::dnodal_position_perturbation_dt_lin_axi_nst(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< DIM >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< DIM >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::PoroelasticityEquations< DIM >::du_dt(), oomph::AdvectionDiffusionEquations< DIM >::du_dt_adv_diff(), oomph::AxisymAdvectionDiffusionEquations::du_dt_axi_adv_diff(), oomph::AxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::AxisymmetricLinearElasticityEquationsBase::du_dt_axisymmetric_linear_elasticity(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::du_dt_cons_adv_diff(), oomph::LinearisedAxisymmetricNavierStokesEquations::du_dt_lin_axi_nst(), oomph::LinearWaveEquations< DIM >::du_dt_lin_wave(), oomph::LinearisedAxisymmetricNavierStokesEquations::du_dt_linearised_axi_nst(), oomph::LinearisedNavierStokesEquations::du_dt_linearised_nst(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::du_dt_nst(), oomph::NavierStokesEquations< DIM >::du_dt_nst(), oomph::PolarNavierStokesEquations::du_dt_pnst(), oomph::SphericalAdvectionDiffusionEquations::du_dt_spherical_adv_diff(), oomph::SphericalNavierStokesEquations::du_dt_spherical_nst(), oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat(), oomph::WomersleyEquations< DIM >::du_dt_womersley(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::extrapolated_strain_rate(), oomph::FiniteElement::get_x(), oomph::ImmersedRigidBodyElement::initialise(), oomph::ImmersedRigidBodyElement::output_centre_of_gravity(), oomph::StraightLine::position(), oomph::Ellipse::position(), oomph::Circle::position(), oomph::PseudoBucklingRing::PseudoBucklingRing(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::SolidICProblem::set_newmark_initial_condition_directly(), oomph::AxisymmetricPoroelasticityEquations::set_q_internal_timestepper(), oomph::PoroelasticityEquations< DIM >::set_q_internal_timestepper(), oomph::SolidICProblem::set_static_initial_condition(), and oomph::PeriodicOrbitEquations::set_timestepper_weights().

◆ time_stepper_pt() [2/2]

TimeStepper* oomph::GeomObject::time_stepper_pt ( ) const
inline

Access function for pointer to time stepper: Null if object is not time-dependent. Const version

200  {
202  }

References Geom_object_time_stepper_pt.

Member Data Documentation

◆ Geom_object_time_stepper_pt

TimeStepper* oomph::GeomObject::Geom_object_time_stepper_pt
protected

◆ Ndim

◆ NLagrangian

unsigned oomph::GeomObject::NLagrangian
protected

Number of Lagrangian (intrinsic) coordinates.

Referenced by nlagrangian(), and set_nlagrangian_and_ndim().


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