oomph::ImmersedRigidBodyTriangleMeshPolygon Class Reference

#include <immersed_rigid_body_elements.h>

+ Inheritance diagram for oomph::ImmersedRigidBodyTriangleMeshPolygon:

Public Member Functions

 ImmersedRigidBodyTriangleMeshPolygon (const Vector< double > &hole_center, const Vector< TriangleMeshCurveSection * > &boundary_polyline_pt, TimeStepper *const &time_stepper_pt, Data *const &centre_displacement_data_pt=0)
 
 ~ImmersedRigidBodyTriangleMeshPolygon ()
 Empty Destuctor. More...
 
void position (const Vector< double > &xi, Vector< double > &r) const
 Overload (again) the position to apply the rotation and translation. More...
 
void position (const unsigned &t, const Vector< double > &xi, Vector< double > &r) const
 Overload (again) the position to apply the rotation and translation. More...
 
void reset_reference_configuration ()
 
- Public Member Functions inherited from oomph::TriangleMeshPolygon
 TriangleMeshPolygon (const Vector< TriangleMeshCurveSection * > &boundary_polyline_pt, const Vector< double > &internal_point_pt=Vector< double >(0), const bool &is_internal_point_fixed=false)
 
virtual ~TriangleMeshPolygon ()
 Empty virtual destructor. More...
 
unsigned ncurve_section () const
 Number of constituent curves. More...
 
unsigned npolyline () const
 Number of constituent polylines. More...
 
TriangleMeshPolyLinepolyline_pt (const unsigned &i) const
 Pointer to i-th constituent polyline. More...
 
TriangleMeshPolyLinepolyline_pt (const unsigned &i)
 Pointer to i-th constituent polyline. More...
 
Vector< unsignedpolygon_boundary_id ()
 Return vector of boundary ids of associated polylines. More...
 
bool is_redistribution_of_segments_between_polylines_enabled ()
 
void enable_redistribution_of_segments_between_polylines ()
 
void disable_redistribution_of_segments_between_polylines ()
 
bool can_update_reference_configuration () const
 Test whether curve can update reference. More...
 
bool is_fixed () const
 Test whether the polygon is fixed or not. More...
 
void set_fixed ()
 Set the polygon to be fixed. More...
 
void set_unfixed ()
 Set the polygon to be allowed to move (default) More...
 
- Public Member Functions inherited from oomph::TriangleMeshClosedCurve
 TriangleMeshClosedCurve (const Vector< TriangleMeshCurveSection * > &curve_section_pt, const Vector< double > &internal_point_pt=Vector< double >(0), const bool &is_internal_point_fixed=false)
 Constructor prototype. More...
 
virtual ~TriangleMeshClosedCurve ()
 Empty destructor. More...
 
unsigned nvertices () const
 Number of vertices. More...
 
unsigned nsegments () const
 Total number of segments. More...
 
void output (std::ostream &outfile, const unsigned &n_sample=50)
 Output each sub-boundary at n_sample (default: 50) points. More...
 
Vector< doubleinternal_point () const
 Coordinates of the internal point. More...
 
Vector< double > & internal_point ()
 Coordinates of the internal point. More...
 
void fix_internal_point ()
 
void unfix_internal_point ()
 
bool is_internal_point_fixed () const
 Test whether the internal point is fixed. More...
 
- Public Member Functions inherited from oomph::TriangleMeshCurve
 TriangleMeshCurve (const Vector< TriangleMeshCurveSection * > &curve_section_pt)
 Empty constructor. More...
 
virtual ~TriangleMeshCurve ()
 Empty destructor. More...
 
unsigned max_boundary_id ()
 Return max boundary id of associated curves. More...
 
void enable_polyline_refinement (const double &tolerance=0.08)
 
void set_polyline_refinement_tolerance (const double &tolerance)
 
double polyline_refinement_tolerance ()
 
void disable_polyline_refinement ()
 Disable refinement of polylines. More...
 
void enable_polyline_unrefinement (const double &tolerance=0.04)
 
void set_polyline_unrefinement_tolerance (const double &tolerance)
 
double polyline_unrefinement_tolerance ()
 
void disable_polyline_unrefinement ()
 Disable unrefinement of polylines. More...
 
virtual TriangleMeshCurveSectioncurve_section_pt (const unsigned &i) const
 Pointer to i-th constituent curve section. More...
 
virtual TriangleMeshCurveSection *& curve_section_pt (const unsigned &i)
 Pointer to i-th constituent curve section. More...
 
- Public Member Functions inherited from oomph::ImmersedRigidBodyElement
 ImmersedRigidBodyElement (GeomObject *const &geom_object_pt, TimeStepper *const &time_stepper_pt, Data *const &centre_displacement_data_pt=0)
 
void set_geometric_rotation ()
 Set the rotation of the object to be included. More...
 
void unset_geometric_rotation ()
 
doubleinitial_phi ()
 Access function for the initial angle. More...
 
doubleinitial_centre_of_mass (const unsigned &i)
 Access function for the initial centre of mass. More...
 
const doubleinitial_centre_of_mass (const unsigned &i) const
 Access function for the initial centre of mass (const version) More...
 
void dposition_dt (const Vector< double > &zeta, const unsigned &j, Vector< double > &drdt)
 Work out the position derivative, including rigid body motion. More...
 
 ~ImmersedRigidBodyElement ()
 Destuctor: Cleanup if required. More...
 
doublemass_shape ()
 
doublemoment_of_inertia_shape ()
 Access to dimensionless polar "moment of inertia" shape parameter. More...
 
Data *& centre_displacement_data_pt ()
 
doublecentre_x_displacement ()
 x-displacement of centre of mass More...
 
doublecentre_y_displacement ()
 y-displacement of centre of mass More...
 
doublecentre_rotation_angle ()
 rotation of centre of mass More...
 
Vector< doublecentre_of_gravity ()
 Get current centre of gravity. More...
 
void pin_centre_of_mass_coordinate (const unsigned &i)
 Pin the i-th coordinate of the centre of mass. More...
 
void unpin_centre_of_mass_coordinate (const unsigned &i)
 Unpin the i-th coordinate of the centre of mass. More...
 
void pin_rotation_angle ()
 Pin the rotation angle. More...
 
void unpin_rotation_angle ()
 Unpin the rotation angle. More...
 
void output_centre_of_gravity (std::ostream &outfile)
 Output position velocity and acceleration of centre of gravity. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Get the contribution to the residuals. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Get residuals including contribution to jacobian. More...
 
void node_update_adjacent_fluid_elements ()
 
void update_in_external_fd (const unsigned &i)
 After an external data change, update the nodal positions. More...
 
void reset_in_external_fd (const unsigned &i)
 Do nothing to reset within finite-differencing of external data. More...
 
void reset_after_external_fd ()
 
void update_in_internal_fd (const unsigned &i)
 After an internal data change, update the nodal positions. More...
 
void reset_in_internal_fd (const unsigned &i)
 Do nothing to reset within finite-differencing of internal data. More...
 
void reset_after_internal_fd ()
 
void get_force_and_torque (const double &time, Vector< double > &force, double &torque)
 
ExternalForceFctPtexternal_force_fct_pt ()
 
ExternalTorqueFctPtexternal_torque_fct_pt ()
 
Mesh *const & drag_mesh_pt ()
 
void set_drag_mesh (Mesh *const &drag_mesh_pt)
 
void flush_drag_mesh ()
 Function to clear the drag mesh and all associated external data. More...
 
unsigned ngeom_data () const
 The position of the object depends on one data item. More...
 
Datageom_data_pt (const unsigned &j)
 
Vector< double > *& g_pt ()
 Access function to the direction of gravity. More...
 
const Vector< double > & g () const
 Access function for gravity. More...
 
double *& re_pt ()
 Access function for the pointer to the fluid Reynolds number. More...
 
const doublere () const
 Access function for the fluid Reynolds number. More...
 
double *& st_pt ()
 Access function for the pointer to the fluid Strouhal number. More...
 
const doublest () const
 Access function for the fluid Strouhal number. More...
 
double *& re_invfr_pt ()
 
const doublere_invfr ()
 Access to the fluid inverse Froude number. More...
 
double *& density_ratio_pt ()
 Access function for the pointer to the density ratio. More...
 
const doubledensity_ratio () const
 Access function for the the density ratio. More...
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual unsigned self_test ()
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- 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 (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
 

Private Member Functions

void get_initial_position (const Vector< double > &xi, Vector< double > &r) const
 Get the initial position of the polygon. More...
 
void assign_zeta ()
 

Private Attributes

Vector< Vector< double > > Zeta_vertex
 Vector of intrisic coordinate values at the nodes. More...
 

Additional Inherited Members

- Public Types inherited from oomph::ImmersedRigidBodyElement
typedef void(* ExternalForceFctPt) (const double &time, Vector< double > &external_force)
 
typedef void(* ExternalTorqueFctPt) (const double &time, double &external_torque)
 
- Static Public Attributes inherited from oomph::GeneralisedElement
static bool Suppress_warning_about_repeated_internal_data
 
static bool Suppress_warning_about_repeated_external_data = true
 
static double Default_fd_jacobian_step = 1.0e-8
 
- Protected Member Functions inherited from oomph::ImmersedRigidBodyElement
 ImmersedRigidBodyElement (TimeStepper *const &time_stepper_pt, Data *const &centre_displacement_data_pt=0)
 
void apply_rigid_body_motion (const unsigned &t, const Vector< double > &initial_x, Vector< double > &r) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void update_before_external_fd ()
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
- Protected Attributes inherited from oomph::TriangleMeshPolygon
bool Enable_redistribution_of_segments_between_polylines
 
bool Can_update_configuration
 
- Protected Attributes inherited from oomph::TriangleMeshClosedCurve
Vector< doubleInternal_point_pt
 Vector of vertex coordinates. More...
 
bool Is_internal_point_fixed
 Indicate whether the internal point should be updated automatically. More...
 
- Protected Attributes inherited from oomph::TriangleMeshCurve
Vector< TriangleMeshCurveSection * > Curve_section_pt
 Vector of curve sections. More...
 
- Protected Attributes inherited from oomph::ImmersedRigidBodyElement
Vector< doubleInitial_centre_of_mass
 X-coordinate of initial centre of gravity. More...
 
double Initial_Phi
 Original rotation angle. More...
 
double Mass
 
double Moment_of_inertia
 Polar moment of inertia of body. More...
 
DataCentre_displacement_data_pt
 
- 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
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

Class upgrading a TriangleMeshPolygon to a "hole" for use during triangle mesh generation. For mesh generation purposes, the main (and only) addition to the base class is the provision of the coordinates of a hole inside the polygon. To faciliate the movement of the "hole" through the domain we also provide a Data object whose three values represent the x and y displacements of its centre of gravity and the polygon's rotation about its centre of gravity. If added to a mesh in the Problem (in its incarnation as a GeneralisedElement) the displacement/rotation of the polygon is computed in response to (i) user-specifiable applied forces and a torque and (ii) the net drag (and associated torque) from a mesh of elements that can exert a drag onto the polygon (typically Navier-Stokes FaceElements that apply a viscous drag to an immersed body, represented by the polygon.)

Constructor & Destructor Documentation

◆ ImmersedRigidBodyTriangleMeshPolygon()

oomph::ImmersedRigidBodyTriangleMeshPolygon::ImmersedRigidBodyTriangleMeshPolygon ( const Vector< double > &  hole_center,
const Vector< TriangleMeshCurveSection * > &  boundary_polyline_pt,
TimeStepper *const &  time_stepper_pt,
Data *const &  centre_displacement_data_pt = 0 
)

Constructor: Specify coordinates of a point inside the hole and a vector of pointers to TriangleMeshPolyLines that define the boundary segments of the polygon. Each TriangleMeshPolyLine has its own boundary ID and can contain multiple (straight-line) segments. The optional final argument is a pointer to a Data object whose three values represent the two displacements of and the rotation angle about the polygon's centre of mass.

392  : TriangleMeshCurve(boundary_polyline_pt),
393  TriangleMeshClosedCurve(boundary_polyline_pt, hole_center),
394  TriangleMeshPolygon(boundary_polyline_pt, hole_center),
396  {
397  // The underlying geometric object can be used to update the configuration
398  // internally before a remesh
399  this->Can_update_configuration = true;
400 
401  // Original rotation angle is zero
402  Initial_Phi = 0.0;
403 
404  // Compute coordinates of centre of gravity etc
405  Vector<double> r_left(2);
406  Vector<double> r_right(2);
407  Mass = 0.0;
408  Initial_centre_of_mass[0] = 0.0;
409  Initial_centre_of_mass[1] = 0.0;
410  double inertia_x = 0.0;
411  double inertia_y = 0.0;
412 
413  // Loop over polylines
414  unsigned nboundary = boundary_polyline_pt.size();
415  for (unsigned i = 0; i < nboundary; i++)
416  {
417  // Loop over the segments to get the vertex coordinates
418  unsigned nseg = boundary_polyline_pt[i]->nsegment();
419  for (unsigned j = 0; j < nseg; j++)
420  {
421  // Get the vertex coordinates
422  r_left = this->polyline_pt(i)->vertex_coordinate(j);
423  r_right = this->polyline_pt(i)->vertex_coordinate(j + 1);
424 
425  // Mass (area)
426  Mass += 0.5 * (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
427 
428  // Centroid
430  (r_left[0] + r_right[0]) *
431  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
433  (r_left[1] + r_right[1]) *
434  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
435  }
436  if (nboundary == 1)
437  {
438  // Get the vertex coordinates
439  r_left = this->polyline_pt(0)->vertex_coordinate(nseg);
440  r_right = this->polyline_pt(0)->vertex_coordinate(0);
441 
442  // Mass (area)
443  Mass += 0.5 * (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
444 
445  // Centroid
447  (r_left[0] + r_right[0]) *
448  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
450  (r_left[1] + r_right[1]) *
451  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
452  }
453  }
454 
455  // Normalise
456  Initial_centre_of_mass[0] /= (6.0 * Mass);
457  Initial_centre_of_mass[1] /= (6.0 * Mass);
458 
459  // Another loop over polylines for moment of inertia
460  for (unsigned i = 0; i < nboundary; i++)
461  {
462  // Loop over the segments to get the vertex coordinates
463  unsigned nseg = boundary_polyline_pt[i]->nsegment();
464  for (unsigned j = 0; j < nseg; j++)
465  {
466  // Get the vertex coordinates
467  r_left = this->polyline_pt(i)->vertex_coordinate(j);
468  r_right = this->polyline_pt(i)->vertex_coordinate(j + 1);
469 
470  // Get moment about centroid
471  r_left[0] -= Initial_centre_of_mass[0];
472  r_left[1] -= Initial_centre_of_mass[1];
473  r_right[0] -= Initial_centre_of_mass[0];
474  r_right[1] -= Initial_centre_of_mass[1];
475 
476  // Moment of inertia
477  inertia_x += 1.0 / 12.0 *
478  (r_left[1] * r_left[1] + r_left[1] * r_right[1] +
479  r_right[1] * r_right[1]) *
480  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
481 
482  inertia_y += 1.0 / 12.0 *
483  (r_left[0] * r_left[0] + r_left[0] * r_right[0] +
484  r_right[0] * r_right[0]) *
485  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
486  }
487 
488  if (nboundary == 1)
489  {
490  // Get the vertex coordinates
491  r_left = this->polyline_pt(0)->vertex_coordinate(nseg);
492  r_right = this->polyline_pt(0)->vertex_coordinate(0);
493 
494  // Get moment about centroid
495  r_left[0] -= Initial_centre_of_mass[0];
496  r_left[1] -= Initial_centre_of_mass[1];
497  r_right[0] -= Initial_centre_of_mass[0];
498  r_right[1] -= Initial_centre_of_mass[1];
499 
500  // Moment of inertia
501  inertia_x += 1.0 / 12.0 *
502  (r_left[1] * r_left[1] + r_left[1] * r_right[1] +
503  r_right[1] * r_right[1]) *
504  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
505 
506  inertia_y += 1.0 / 12.0 *
507  (r_left[0] * r_left[0] + r_left[0] * r_right[0] +
508  r_right[0] * r_right[0]) *
509  (r_left[0] * r_right[1] - r_right[0] * r_left[1]);
510  }
511  }
512 
513  // Polar moment of inertia is sum of two orthogonal planar moments
514  Moment_of_inertia = inertia_x + inertia_y;
515 
516  // // Tested for circular and elliptical cross section
517  // cout << "Mass : " << Mass << std::endl;
518  // cout << "Moment of inertia: " << Moment_of_inertia << std::endl;
519  // cout << "X_c : " << Initial_centre_of_mass[0] <<
520  // std::endl; cout << "Y_c : " << Initial_centre_of_mass[1]
521  // << std::endl; pause("done");
522 
523 
524  // Assign the intrinsic coordinate
525  this->assign_zeta();
526 
527  // {
528  // unsigned n_poly = this->npolyline();
529  // for(unsigned p=0;p<n_poly;++p)
530  // {
531  // std::cout << "Polyline " << p << "\n";
532  // std::cout << "-----------------------\n";
533  // unsigned n_vertex = Zeta_vertex[p].size();
534  // for(unsigned v=0;v<n_vertex;v++)
535  // {
536  // std::cout << v << " " << Zeta_vertex[p][v] << "\n";
537  // }
538  // }
539  // }
540  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
ImmersedRigidBodyElement(TimeStepper *const &time_stepper_pt, Data *const &centre_displacement_data_pt=0)
Definition: immersed_rigid_body_elements.h:110
double Initial_Phi
Original rotation angle.
Definition: immersed_rigid_body_elements.h:582
double Moment_of_inertia
Polar moment of inertia of body.
Definition: immersed_rigid_body_elements.h:588
Data *& centre_displacement_data_pt()
Definition: immersed_rigid_body_elements.h:235
double Mass
Definition: immersed_rigid_body_elements.h:585
Vector< double > Initial_centre_of_mass
X-coordinate of initial centre of gravity.
Definition: immersed_rigid_body_elements.h:579
void assign_zeta()
Definition: immersed_rigid_body_elements.h:806
TriangleMeshClosedCurve(const Vector< TriangleMeshCurveSection * > &curve_section_pt, const Vector< double > &internal_point_pt=Vector< double >(0), const bool &is_internal_point_fixed=false)
Constructor prototype.
Definition: unstructured_two_d_mesh_geometry_base.cc:1420
TriangleMeshCurve(const Vector< TriangleMeshCurveSection * > &curve_section_pt)
Empty constructor.
Definition: unstructured_two_d_mesh_geometry_base.h:1139
Vector< double > vertex_coordinate(const unsigned &i) const
Coordinate vector of i-th vertex (const version)
Definition: unstructured_two_d_mesh_geometry_base.h:929
TriangleMeshPolyLine * polyline_pt(const unsigned &i) const
Pointer to i-th constituent polyline.
Definition: unstructured_two_d_mesh_geometry_base.h:1488
TriangleMeshPolygon(const Vector< TriangleMeshCurveSection * > &boundary_polyline_pt, const Vector< double > &internal_point_pt=Vector< double >(0), const bool &is_internal_point_fixed=false)
Definition: unstructured_two_d_mesh_geometry_base.cc:1577
bool Can_update_configuration
Definition: unstructured_two_d_mesh_geometry_base.h:1624
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References assign_zeta(), oomph::TriangleMeshPolygon::Can_update_configuration, i, oomph::ImmersedRigidBodyElement::Initial_centre_of_mass, oomph::ImmersedRigidBodyElement::Initial_Phi, j, oomph::ImmersedRigidBodyElement::Mass, oomph::ImmersedRigidBodyElement::Moment_of_inertia, oomph::TriangleMeshPolygon::polyline_pt(), and oomph::TriangleMeshPolyLine::vertex_coordinate().

◆ ~ImmersedRigidBodyTriangleMeshPolygon()

oomph::ImmersedRigidBodyTriangleMeshPolygon::~ImmersedRigidBodyTriangleMeshPolygon ( )
inline

Empty Destuctor.

681 {}

Member Function Documentation

◆ assign_zeta()

void oomph::ImmersedRigidBodyTriangleMeshPolygon::assign_zeta ( )
inlineprivate

Helper function to assign the values of the (scaled) arc-length to each node of each polyline. The direction will be the natural order of the vertices within the polyline.

807  {
808  // Find the number of polylines
809  unsigned n_poly = this->npolyline();
810 
811  // Allocate appropriate storage for the zeta values
812  Zeta_vertex.resize(n_poly);
813 
814  // Temporary storage for the vertex coordinates
815  Vector<double> vertex_coord_first;
816  Vector<double> vertex_coord_next;
817 
818  // Set the initial value of zeta
819  double zeta_offset = 0.0;
820 
821  // Loop over the polylines
822  for (unsigned p = 0; p < n_poly; ++p)
823  {
824  // Cache the pointer to the polyline
825  TriangleMeshPolyLine* const line_pt = this->polyline_pt(p);
826 
827  // Find the number of vertices in the polyline
828  unsigned n_vertex = line_pt->nvertex();
829 
830  // Allocate storage and set initial value
831  Zeta_vertex[p].resize(n_vertex);
832  Zeta_vertex[p][0] = 0.0;
833 
834  // Loop over the vertices in the polyline and calculate the length
835  // between each for use as the intrinsic coordinate
836  vertex_coord_first = line_pt->vertex_coordinate(0);
837  for (unsigned v = 1; v < n_vertex; v++)
838  {
839  vertex_coord_next = line_pt->vertex_coordinate(v);
840  double length =
841  sqrt(pow(vertex_coord_next[0] - vertex_coord_first[0], 2.0) +
842  pow(vertex_coord_next[1] - vertex_coord_first[1], 2.0));
843  Zeta_vertex[p][v] = Zeta_vertex[p][v - 1] + length;
844  vertex_coord_first = vertex_coord_next;
845  }
846 
847  // Now scale the length to unity and add the offset
848  Zeta_vertex[p][0] += zeta_offset;
849  for (unsigned v = 1; v < n_vertex; v++)
850  {
851  Zeta_vertex[p][v] /= Zeta_vertex[p][n_vertex - 1];
852  Zeta_vertex[p][v] += zeta_offset;
853  }
854  zeta_offset += 1.0;
855  } // End of loop over polylines
856  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
float * p
Definition: Tutorial_Map_using.cpp:9
Vector< Vector< double > > Zeta_vertex
Vector of intrisic coordinate values at the nodes.
Definition: immersed_rigid_body_elements.h:859
unsigned npolyline() const
Number of constituent polylines.
Definition: unstructured_two_d_mesh_geometry_base.h:1482
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625

References oomph::TriangleMeshPolygon::npolyline(), oomph::TriangleMeshPolyLine::nvertex(), p, oomph::TriangleMeshPolygon::polyline_pt(), Eigen::bfloat16_impl::pow(), sqrt(), v, oomph::TriangleMeshPolyLine::vertex_coordinate(), and Zeta_vertex.

Referenced by ImmersedRigidBodyTriangleMeshPolygon().

◆ get_initial_position()

void oomph::ImmersedRigidBodyTriangleMeshPolygon::get_initial_position ( const Vector< double > &  xi,
Vector< double > &  r 
) const
inlineprivate

Get the initial position of the polygon.

712  {
713  // Find the number of polylines (boundaries)
714  unsigned n_poly = this->npolyline();
715 
716  // The boundary coordinate will be contiguous from polyline to
717  // polyline and each polyline will have the scaled arclength coordinate
718  // in the range 0->1.
719 
720  // Find the maximum coordinate
721  double zeta_max = Zeta_vertex[n_poly - 1].back();
722 
723  // If we are above the maximum complain
724  if (xi[0] > zeta_max)
725  {
726  std::ostringstream error_message;
727  error_message << "Value of intrinsic coordinate " << xi[0]
728  << "greater than maximum " << zeta_max << "\n";
729  throw OomphLibError(error_message.str(),
730  "TriangleMeshPolygon::position()",
732  }
733 
734  // If equal to the maximum, return the final vertex
735  if (xi[0] == zeta_max)
736  {
737  unsigned n_vertex = this->polyline_pt(n_poly - 1)->nvertex();
738  r = this->polyline_pt(n_poly - 1)->vertex_coordinate(n_vertex - 1);
739  return;
740  }
741 
742  // Otherwise
743 
744  // Find out which polyline we are in
745  // If we've got here this should be less than n_poly
746  unsigned p = static_cast<unsigned>(floor(xi[0]));
747 
748  // If we are "above" the last polyline then throw an error
749  // This should have been caught by the trap above
750  if (p >= n_poly)
751  {
752  std::ostringstream error_message;
753  error_message
754  << "Something has gone wrong.\n"
755  << "The integer part of the input intrinsic coordinate is " << p
756  << "\nwhich is equal to or greater than the number of polylines, "
757  << n_poly << ".\n"
758  << "This should have triggered an earlier error\n";
759 
760 
761  throw OomphLibError(error_message.str(),
764  }
765 
766  // Cache the appropriate polyline
767  TriangleMeshPolyLine* const line_pt = this->polyline_pt(p);
768 
769  // If we are at the first vertex in the polyline, return it
770  if (xi[0] == Zeta_vertex[p][0])
771  {
772  r = line_pt->vertex_coordinate(0);
773  return;
774  }
775 
776  // Otherwise loop over the other points to find the appropriate
777  // segment
778 
779  // Find the number of vertices in the chosen polyline
780  unsigned n_vertex = line_pt->nvertex();
781  // Now start from the first node in the appropriate polyline and loop up
782  for (unsigned v = 1; v < n_vertex; v++)
783  {
784  // First time that zeta falls below the vertex coordinate
785  // we have something
786  if (xi[0] < Zeta_vertex[p][v])
787  {
788  double fraction = (xi[0] - Zeta_vertex[p][v - 1]) /
789  (Zeta_vertex[p][v] - Zeta_vertex[p][v - 1]);
790  Vector<double> first = line_pt->vertex_coordinate(v - 1);
791  Vector<double> last = line_pt->vertex_coordinate(v);
792  r.resize(2);
793  for (unsigned i = 0; i < 2; i++)
794  {
795  r[i] = first[i] + fraction * (last[i] - first[i]);
796  }
797  return;
798  }
799  }
800  }
unsigned nvertex() const
Number of vertices.
Definition: unstructured_two_d_mesh_geometry_base.h:904
static constexpr const last_t last
Definition: IndexedViewHelper.h:48
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 floor(const bfloat16 &a)
Definition: BFloat16.h:643
r
Definition: UniformPSDSelfTest.py:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Eigen::bfloat16_impl::floor(), i, Eigen::placeholders::last, oomph::TriangleMeshPolygon::npolyline(), oomph::TriangleMeshPolyLine::nvertex(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, p, oomph::TriangleMeshPolygon::polyline_pt(), UniformPSDSelfTest::r, v, oomph::TriangleMeshPolyLine::vertex_coordinate(), and Zeta_vertex.

Referenced by position().

◆ position() [1/2]

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

Overload (again) the position to apply the rotation and translation.

Reimplemented from oomph::ImmersedRigidBodyElement.

696  {
697  Vector<double> initial_x(2);
698  this->get_initial_position(xi, initial_x);
699  this->apply_rigid_body_motion(t, initial_x, r);
700  }
void apply_rigid_body_motion(const unsigned &t, const Vector< double > &initial_x, Vector< double > &r) const
Definition: immersed_rigid_body_elements.h:498
void get_initial_position(const Vector< double > &xi, Vector< double > &r) const
Get the initial position of the polygon.
Definition: immersed_rigid_body_elements.h:711

References oomph::ImmersedRigidBodyElement::apply_rigid_body_motion(), get_initial_position(), and UniformPSDSelfTest::r.

◆ position() [2/2]

void oomph::ImmersedRigidBodyTriangleMeshPolygon::position ( const Vector< double > &  xi,
Vector< double > &  r 
) const
inlinevirtual

Overload (again) the position to apply the rotation and translation.

Reimplemented from oomph::ImmersedRigidBodyElement.

685  {
686  Vector<double> initial_x(2);
687  this->get_initial_position(xi, initial_x);
688  this->apply_rigid_body_motion(0, initial_x, r);
689  }

References oomph::ImmersedRigidBodyElement::apply_rigid_body_motion(), get_initial_position(), and UniformPSDSelfTest::r.

◆ reset_reference_configuration()

void oomph::ImmersedRigidBodyTriangleMeshPolygon::reset_reference_configuration ( )
virtual

Update the reference configuration by re-setting the original position of the vertices to their current ones, re-set the original position of the centre of mass, and the displacements and rotations relative to it

Reimplemented from oomph::TriangleMeshPolygon.

550  {
551  Vector<double> x_orig(2);
552  Vector<double> r(2);
553 
554  // Loop over the polylines and update their vertex positions
555  unsigned npoly = this->ncurve_section();
556  for (unsigned i = 0; i < npoly; i++)
557  {
558  TriangleMeshPolyLine* poly_line_pt = this->polyline_pt(i);
559  unsigned nvertex = poly_line_pt->nvertex();
560  for (unsigned j = 0; j < nvertex; j++)
561  {
562  x_orig = poly_line_pt->vertex_coordinate(j);
563  this->apply_rigid_body_motion(0, x_orig, r);
564  poly_line_pt->vertex_coordinate(j) = r;
565  }
566  }
567 
568  // Update coordinates of hole
569  Vector<double> orig_hole_coord(this->internal_point());
570  this->apply_rigid_body_motion(0, orig_hole_coord, this->internal_point());
571 
572  // Update centre of gravity
573  double x_displ = Centre_displacement_data_pt->value(0);
574  double y_displ = Centre_displacement_data_pt->value(1);
575  double phi_displ = Centre_displacement_data_pt->value(2);
576  Initial_centre_of_mass[0] += x_displ;
577  Initial_centre_of_mass[1] += y_displ;
578  Initial_Phi += phi_displ;
579 
580  // Reset displacement and rotation ("non-previous-value"
581  // history values stay)
582  TimeStepper* timestepper_pt =
584  unsigned nprev = timestepper_pt->nprev_values();
585  for (unsigned t = 0; t < nprev; t++)
586  {
588  t, 0, Centre_displacement_data_pt->value(t, 0) - x_displ);
589 
591  t, 1, Centre_displacement_data_pt->value(t, 1) - y_displ);
592 
594  t, 2, Centre_displacement_data_pt->value(t, 2) - phi_displ);
595  }
596  }
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
double value(const unsigned &i) const
Definition: nodes.h:293
Data * Centre_displacement_data_pt
Definition: immersed_rigid_body_elements.h:592
virtual unsigned nprev_values() const =0
Number of previous values available: 0 for static, 1 for BDF<1>,...
Vector< double > internal_point() const
Coordinates of the internal point.
Definition: unstructured_two_d_mesh_geometry_base.h:1403
unsigned ncurve_section() const
Number of constituent curves.
Definition: unstructured_two_d_mesh_geometry_base.h:1476
t
Definition: plotPSD.py:36

References oomph::ImmersedRigidBodyElement::apply_rigid_body_motion(), oomph::ImmersedRigidBodyElement::Centre_displacement_data_pt, i, oomph::ImmersedRigidBodyElement::Initial_centre_of_mass, oomph::ImmersedRigidBodyElement::Initial_Phi, oomph::TriangleMeshClosedCurve::internal_point(), j, oomph::TriangleMeshPolygon::ncurve_section(), oomph::TimeStepper::nprev_values(), oomph::TriangleMeshPolyLine::nvertex(), oomph::TriangleMeshPolygon::polyline_pt(), UniformPSDSelfTest::r, oomph::Data::set_value(), plotPSD::t, oomph::Data::time_stepper_pt(), oomph::Data::value(), and oomph::TriangleMeshPolyLine::vertex_coordinate().

Member Data Documentation

◆ Zeta_vertex

Vector<Vector<double> > oomph::ImmersedRigidBodyTriangleMeshPolygon::Zeta_vertex
private

Vector of intrisic coordinate values at the nodes.

Referenced by assign_zeta(), and get_initial_position().


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