CircularPenetratorElement Class Reference
+ Inheritance diagram for CircularPenetratorElement:

Public Member Functions

 CircularPenetratorElement (SolidNode *control_node_pt, const unsigned &index_of_contact_pressure, const unsigned &index_of_horizontal_displacement, const unsigned &index_of_vertical_displacement, double *r_pt)
 
Vector< std::pair< Data *, unsigned > > equilibrium_data ()
 
double angle () const
 Angle of rotation around contact point. More...
 
void set_angle (const double &angle)
 Set angle of rotation around contact point. More...
 
Meshcontact_element_mesh_pt () const
 
void set_contact_element_mesh_pt (Mesh *contact_element_mesh_pt)
 
void set_equilibrium_target_forces ()
 Set target horizontal and vertical force to be in equilibrium. More...
 
double target_weight ()
 Target weight (returns zero if not imposed) More...
 
double target_horizontal_force ()
 Target horizontal force (returns zero if not imposed) More...
 
double target_yc ()
 Target vertical position of control point (returns zero if not imposed) More...
 
double target_rotation_angle ()
 Target rotation angle about contact point (returns zero if not imposed) More...
 
bool yc_is_imposed ()
 Is vertical positon of control node imposed? If false then weight imposed. More...
 
void impose_weight (double *target_weight_pt)
 
void impose_yc (double *target_yc_pt)
 
bool rotation_angle_is_imposed ()
 
void impose_horizontal_force (double *target_horizontal_force_pt)
 
void impose_rotation_angle (double *target_rotation_angle_pt)
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Fill in contribution to residuals. More...
 
Vector< doublecentre () const
 Get centre of penetrator. More...
 
double centre (const unsigned &i) const
 Get centre of penetrator. More...
 
void penetration (const Vector< double > &x, const Vector< double > &n, double &d, bool &intersection) const
 Get penetration for given point x. More...
 
void output (std::ostream &outfile, const unsigned &nplot) const
 Output coordinates of penetrator at nplot plot points. More...
 
Vector< doubleresulting_force () const
 Resulting force from all associated ContactElements. More...
 
double radius () const
 Radius of penetrator. More...
 
 CircularPenetratorElement (SolidNode *control_node_pt, const unsigned &index_of_contact_pressure, double *r_pt)
 
virtual ~CircularPenetratorElement ()
 Virtual destructor. More...
 
Vector< std::pair< Data *, unsigned > > equilibrium_data ()
 
double angle () const
 Angle of rotation around contact point. More...
 
void set_angle (const double &angle)
 Set angle of rotation around contact point. More...
 
Meshcontact_element_mesh_pt () const
 
void set_contact_element_mesh_pt (Mesh *contact_element_mesh_pt)
 
void set_equilibrium_target_forces ()
 Set target horizontal and vertical force to be in equilibrium. More...
 
double target_weight ()
 Target weight (returns zero if not imposed) More...
 
double target_horizontal_force ()
 Target horizontal force (returns zero if not imposed) More...
 
double target_yc ()
 Target vertical position of control point (returns zero if not imposed) More...
 
double target_rotation_angle ()
 Target rotation angle about contact point (returns zero if not imposed) More...
 
bool yc_is_imposed ()
 Is vertical positon of control node imposed? If false then weight imposed. More...
 
void impose_weight (double *target_weight_pt)
 
void impose_yc (double *target_yc_pt)
 
bool rotation_angle_is_imposed ()
 
void impose_horizontal_force (double *target_horizontal_force_pt)
 
void impose_rotation_angle (double *target_rotation_angle_pt)
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Fill in contribution to residuals. More...
 
Vector< doublecentre () const
 Get centre of penetrator. More...
 
double centre (const unsigned &i) const
 Get centre of penetrator. More...
 
void penetration (const Vector< double > &x, const Vector< double > &n, double &d, bool &intersection) const
 Get penetration for given point x. More...
 
void output (std::ostream &outfile, const unsigned &nplot) const
 Output coordinates of penetrator at nplot plot points. More...
 
Vector< doubleresulting_force () const
 Resulting force from all associated ContactElements. More...
 
double radius () const
 Radius of penetrator. 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::Penetrator
 Penetrator ()
 Constructor. More...
 
virtual ~Penetrator ()
 Destructor. More...
 
virtual Vector< doublerigid_body_displacement () const
 
virtual void surface_coordinate (const Vector< double > &x, Vector< double > &zeta) const
 

Private Attributes

doubleRadius_pt
 Pointer to radius of penetrator. More...
 
SolidNodeControl_node_pt
 Control node. More...
 
unsigned Index_of_contact_pressure
 
unsigned Index_of_vertical_displacement
 Where is the vertical displacement (linear_elasticity) stored? More...
 
unsigned Index_of_horizontal_displacement
 Where is the horizontal displacement (linear_elasticity) stored? More...
 
unsigned External_data_index_of_traded_contact_pressure
 
doubleTarget_weight_pt
 
doubleTarget_horizontal_force_pt
 
doubleTarget_yc_pt
 
doubleTarget_rotation_angle_pt
 
MeshContact_element_mesh_pt
 Mesh of contact elements that contribute to weight/horizontal force. More...
 

Additional Inherited Members

- 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::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 reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
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)
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////
Penetrator that keeps circle in contact with a control node on target surface (made of solid contact face elements) – centre of the circular penetrator is located at

{\bf r}_c = {\bf R}_p + R {\bf e}_alpha

where {\bf R}_p is the position of the control point, R the radius of the circular penetrator, and {\bf e}_alpha is a unit vector inclined at an angle \alpha against the vertical. Penetration can be driven in two ways. (1) We impose the vertical position of the control point (by pseudo-hijacking the Lagrange-multiplier (representing the contact pressure) stored at the controlled node. This means that rather than determining the contact pressure from the no-penetration constraint, (which we know to be satisfied by construction) we determine it from the condition that {\bf R}_p \cdot {\bf e}_y = Y_c which is prescribed. We also impose the angle \alpha (stored as an internal Data value in the element) by solving it via the equation \alpha-\alpha_{prescribed} = 0. (2) We impose the weight (i.e. the vertical reaction force from the contact elements) by using the equation

\int p_c {\bf n} \cdot {\bf e}_y ds - W = 0

as the equation for the pseudo-hijacked contact pressure and similarly, use the horizontal force balance

\int p_c {\bf n} \cdot {\bf e}_x ds - H = 0

to determine the rotation angle. Here, W and H are prescribed and the integral is computed from the contact elements that potentially impact on the penetrator.

Constructor & Destructor Documentation

◆ CircularPenetratorElement() [1/2]

CircularPenetratorElement::CircularPenetratorElement ( SolidNode control_node_pt,
const unsigned index_of_contact_pressure,
const unsigned index_of_horizontal_displacement,
const unsigned index_of_vertical_displacement,
double r_pt 
)
inline

Constructor: Pass pointer to control node whose index_of_contact_pressure-th value represents the Lagrange multiplier (the discrete contact pressure) that has been traded for the vertical displacement/weight constraint. Also need the indices of the nodal values that store the horizontal/vertical displacement (linear elasticity).

272  {
273  // Create internal data, representing the angle of rotation about
274  // contact point. Determined either directly via insisting
275  // that the difference between this value and the target is zero
276  // or by insisting that the horizontal force is zero.
277  add_internal_data(new Data(1));
278  internal_data_pt(0)->set_value(0,0.0);
279 
280  // Store pointer to radius
281  Radius_pt=r_pt;
282 
283  // Control node
284  Control_node_pt=control_node_pt;
285 
286  // Where is the tradedd contact pressure stored?
287  Index_of_contact_pressure=index_of_contact_pressure;
288 
289  // Where is the horizontal displacement (linear_elasticity) stored?
290  Index_of_horizontal_displacement=index_of_horizontal_displacement;
291 
292  // Where is the vertical displacement (linear_elasticity) stored?
293  Index_of_vertical_displacement=index_of_vertical_displacement;
294 
295  //Pointer to target weight (null if vertical displacement of control
296  //node is imposed)
298 
299  // Pointer to target horizontal force (null if rotation angle angle
300  // about control node is imposed)
302 
303  // Pointer to target vertical displacement of control node (null if
304  // weight is imposed)
305  Target_yc_pt=0;
306 
307  // Pointer to target rotation angle about control node (null
308  // if horizontal force is imposed)
310 
311  // Pointer to mesh of contact elements that contribute to force
313  }
unsigned Index_of_horizontal_displacement
Where is the horizontal displacement (linear_elasticity) stored?
Definition: linear_solid_contact_with_gravity.cc:769
unsigned Index_of_contact_pressure
Definition: linear_solid_contact_with_gravity.cc:763
unsigned Index_of_vertical_displacement
Where is the vertical displacement (linear_elasticity) stored?
Definition: linear_solid_contact_with_gravity.cc:766
double * Target_horizontal_force_pt
Definition: linear_solid_contact_with_gravity.cc:781
double * Target_yc_pt
Definition: linear_solid_contact_with_gravity.cc:785
SolidNode * Control_node_pt
Control node.
Definition: linear_solid_contact_with_gravity.cc:759
double * Target_weight_pt
Definition: linear_solid_contact_with_gravity.cc:777
double * Target_rotation_angle_pt
Definition: linear_solid_contact_with_gravity.cc:789
Mesh * Contact_element_mesh_pt
Mesh of contact elements that contribute to weight/horizontal force.
Definition: linear_solid_contact_with_gravity.cc:792
double * Radius_pt
Pointer to radius of penetrator.
Definition: linear_solid_contact_with_gravity.cc:756
Definition: nodes.h:86
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:62

◆ CircularPenetratorElement() [2/2]

CircularPenetratorElement::CircularPenetratorElement ( SolidNode control_node_pt,
const unsigned index_of_contact_pressure,
double r_pt 
)
inline

Constructor: Pass pointer to control node whose index_of_contact_pressure-th value represents the Lagrange multiplier (the discrete contact pressure) that has been traded for the vertical displacement/weight constraint.

106  {
107  // Create internal data, representing the angle of rotation about
108  // contact point. Determined either directly via insisting
109  // that the difference between this value and the target is zero
110  // or by insisting that the horizontal force is zero.
111  add_internal_data(new Data(1));
112  internal_data_pt(0)->set_value(0,0.0);
113 
114  // Store pointer to radius
115  Radius_pt=r_pt;
116 
117  // Control node
118  Control_node_pt=control_node_pt;
119 
120  // Where is the tradedd contact pressure stored?
121  Index_of_contact_pressure=index_of_contact_pressure;
122 
123  //Pointer to target weight (null if vertical displacement of control
124  //node is imposed)
126 
127  // Pointer to target horizontal force (null if rotation angle angle
128  // about control node is imposed)
130 
131  // Pointer to target vertical displacement of control node (null if
132  // weight is imposed)
133  Target_yc_pt=0;
134 
135  // Pointer to target rotation angle about control node (null
136  // if horizontal force is imposed)
138 
139  // Pointer to mesh of contact elements that contribute to force
141  }

◆ ~CircularPenetratorElement()

virtual CircularPenetratorElement::~CircularPenetratorElement ( )
inlinevirtual

Virtual destructor.

145 {};

Member Function Documentation

◆ angle() [1/2]

double CircularPenetratorElement::angle ( ) const
inline

Angle of rotation around contact point.

355  {
356  return internal_data_pt(0)->value(0);
357  }
double value(const unsigned &i) const
Definition: nodes.h:293

Referenced by ContactProblem< ELEMENT >::complete_problem_setup().

◆ angle() [2/2]

double CircularPenetratorElement::angle ( ) const
inline

Angle of rotation around contact point.

186  {
187  return internal_data_pt(0)->value(0);
188  }

◆ centre() [1/4]

Vector<double> CircularPenetratorElement::centre ( ) const
inline

Get centre of penetrator.

647  {
648  Vector<double> rc(2);
649  rc[0]=centre(0);
650  rc[1]=centre(1);
651  return rc;
652  }
Vector< double > centre() const
Get centre of penetrator.
Definition: linear_solid_contact_with_gravity.cc:646
list rc
Definition: plotDoE.py:16

References plotDoE::rc.

◆ centre() [2/4]

Vector<double> CircularPenetratorElement::centre ( ) const
inline

Get centre of penetrator.

470  {
471  Vector<double> rc(2);
472  rc[0]=centre(0);
473  rc[1]=centre(1);
474  return rc;
475  }

References plotDoE::rc.

◆ centre() [3/4]

double CircularPenetratorElement::centre ( const unsigned i) const
inline

Get centre of penetrator.

656  {
657  switch (i)
658  {
659  case 0:
660  return Control_node_pt->x(0)+
662  (*Radius_pt)*sin(angle());
663  break;
664 
665  case 1:
666  return Control_node_pt->x(1)+
668  (*Radius_pt)*cos(angle());
669  break;
670 
671  default:
672  std::stringstream junk;
673  junk << "Wrong index: " << i
674  << "\nCan only handle 0 or 1 (it's 2D!)\n";
675  throw OomphLibError(
676  junk.str(),
679  }
680  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137
int i
Definition: BiCGSTAB_step_by_step.cpp:9
double angle() const
Angle of rotation around contact point.
Definition: linear_solid_contact_with_gravity.cc:354
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
double value(const unsigned &i) const
Definition: nodes.cc:2408
Definition: oomph_definitions.h:222
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Jeffery_Solution::angle(), cos(), i, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and sin().

◆ centre() [4/4]

double CircularPenetratorElement::centre ( const unsigned i) const
inline

Get centre of penetrator.

479  {
480  switch (i)
481  {
482  case 0:
483  return Control_node_pt->x(0)+(*Radius_pt)*sin(angle());
484  break;
485 
486  case 1:
487  return Control_node_pt->x(1)+(*Radius_pt)*cos(angle());
488  break;
489 
490  default:
491  std::stringstream junk;
492  junk << "Wrong index: " << i
493  << "\nCan only handle 0 or 1 (it's 2D!)\n";
494  throw OomphLibError(
495  junk.str(),
498  }
499  }

References Jeffery_Solution::angle(), cos(), i, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and sin().

◆ contact_element_mesh_pt() [1/2]

Mesh* CircularPenetratorElement::contact_element_mesh_pt ( ) const
inline

Access to pointer to mesh of contact elements that contribute to force on penetrator

370  {
372  }

◆ contact_element_mesh_pt() [2/2]

Mesh* CircularPenetratorElement::contact_element_mesh_pt ( ) const
inline

Access to pointer to mesh of contact elements that contribute to force on penetrator

201  {
203  }

◆ equilibrium_data() [1/2]

Vector<std::pair<Data*,unsigned> > CircularPenetratorElement::equilibrium_data ( )
inlinevirtual

Vector of pairs identifying values (via a pair of pointer to Data object and index within it) that correspond to the Data values that are determined by the horizontal/vertical/... equilibrium equations.

Reimplemented from oomph::Penetrator.

320  {
321  // We're in 2D
323 
324 
325  // Horizontal equilibrium determines the rotation angle
326  // which is stored as the zero-th internal data
328  {
329  thingy[0]=std::make_pair(static_cast<Data*>(0),0);
330  }
331  else
332  {
333  thingy[0]=std::make_pair(internal_data_pt(0),0);
334  }
335 
336 
337  // Vertical equilibrium determines the discrete contact pressure
338  // (Lagrange multiplier) at control node
339  if (Target_weight_pt==0)
340  {
341  thingy[1]=std::make_pair(static_cast<Data*>(0),0);
342  }
343  else
344  {
345  thingy[1]=std::make_pair(
348  }
349 
350  return thingy;
351  }
unsigned External_data_index_of_traded_contact_pressure
Definition: linear_solid_contact_with_gravity.cc:773
Data *& external_data_pt(const unsigned &i)
Return a pointer to i-th external data object.
Definition: elements.h:659

◆ equilibrium_data() [2/2]

Vector<std::pair<Data*,unsigned> > CircularPenetratorElement::equilibrium_data ( )
inlinevirtual

Vector of pairs identifying values (via a pair of pointer to Data object and index within it) that correspond to the Data values that are determined by the horizontal/vertical/... equilibrium equations.

Reimplemented from oomph::Penetrator.

151  {
152  // We're in 2D
154 
155 
156  // Horizontal equilibrium determines the rotation angle
157  // which is stored as the zero-th internal data
159  {
160  thingy[0]=std::make_pair(static_cast<Data*>(0),0);
161  }
162  else
163  {
164  thingy[0]=std::make_pair(internal_data_pt(0),0);
165  }
166 
167 
168  // Vertical equilibrium determines the discrete contact pressure
169  // (Lagrange multiplier) at control node
170  if (Target_weight_pt==0)
171  {
172  thingy[1]=std::make_pair(static_cast<Data*>(0),0);
173  }
174  else
175  {
176  thingy[1]=std::make_pair(
179  }
180 
181  return thingy;
182  }

◆ fill_in_contribution_to_residuals() [1/2]

void CircularPenetratorElement::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Fill in contribution to residuals.

Reimplemented from oomph::GeneralisedElement.

565  {
566  // Get resulting force from all associated PseudoContactElements
567  // onto the elastic body
569 
570  // Equation for Lagrange multiplier (contact pressure) at controlled
571  // node
572  int local_eqn=external_local_eqn(
574  if (local_eqn>=0)
575  {
576  // Resulting force from all associated PseudoContactElements
577  // onto the elastic body is equal and opposite to force on penetrator
578  if (Target_weight_pt!=0)
579  {
580  residuals[local_eqn]+=(*Target_weight_pt);
581  }
582  // Impose vertical position of control node
583  else
584  {
585 #ifdef PARANOID
586  if (Target_yc_pt!=0)
587  {
588 #endif
589  residuals[local_eqn]+=Control_node_pt->x(1)+
591 #ifdef PARANOID
592  }
593  else
594  {
595  std::stringstream junk;
596  junk << "Target_yc_pt=0\n"
597  << "Set with impose_yc(...)\n";
598  throw OomphLibError(
599  junk.str(),
602  }
603 #endif
604  }
605  }
606 
607 
608 
609  // Equation for rotation angle
610  local_eqn=internal_local_eqn(0,0);
611  if (local_eqn>=0)
612  {
613  // Resulting force from all associated PseudoContactElements
614  // onto the elastic body is equal and opposite to force on penetrator
616  {
617  residuals[local_eqn]+=(*Target_horizontal_force_pt);
618  }
619  // Set rotation angle
620  else
621  {
622 #ifdef PARANOID
624  {
625 #endif
626  residuals[local_eqn]+=internal_data_pt(0)->value(0)-
627  (*Target_rotation_angle_pt);
628 #ifdef PARANOID
629  }
630  else
631  {
632  std::stringstream junk;
633  junk << "Target_rotation_angle_pt=0\n"
634  << "Set with impose_rotation_angle(...)\n";
635  throw OomphLibError(
636  junk.str(),
639  }
640 #endif
641  }
642  }
643  }
Vector< double > resulting_force() const
Resulting force from all associated ContactElements.
Definition: linear_solid_contact_with_gravity.cc:732
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Definition: elements.h:267
int external_local_eqn(const unsigned &i, const unsigned &j)
Definition: elements.h:311

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ fill_in_contribution_to_residuals() [2/2]

void CircularPenetratorElement::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Fill in contribution to residuals.

Reimplemented from oomph::GeneralisedElement.

389  {
390  // Get resulting force from all associated PseudoContactElements
391  // onto the elastic body
393 
394  // Equation for Lagrange multiplier (contact pressure) at controlled
395  // node
396  int local_eqn=external_local_eqn(
398  if (local_eqn>=0)
399  {
400  // Resulting force from all associated PseudoContactElements
401  // onto the elastic body is equal and opposite to force on penetrator
402  if (Target_weight_pt!=0)
403  {
404  residuals[local_eqn]+=(*Target_weight_pt);
405  }
406  // Impose vertical position of control node
407  else
408  {
409 #ifdef PARANOID
410  if (Target_yc_pt!=0)
411  {
412 #endif
413  residuals[local_eqn]+=Control_node_pt->x(1)-(*Target_yc_pt);
414 #ifdef PARANOID
415  }
416  else
417  {
418  std::stringstream junk;
419  junk << "Target_yc_pt=0\n"
420  << "Set with impose_yc(...)\n";
421  throw OomphLibError(
422  junk.str(),
425  }
426 #endif
427  }
428  }
429 
430 
431 
432  // Equation for rotation angle
433  local_eqn=internal_local_eqn(0,0);
434  if (local_eqn>=0)
435  {
436  // Resulting force from all associated PseudoContactElements
437  // onto the elastic body is equal and opposite to force on penetrator
439  {
440  residuals[local_eqn]+=(*Target_horizontal_force_pt);
441  }
442  // Set rotation angle
443  else
444  {
445 #ifdef PARANOID
447  {
448 #endif
449  residuals[local_eqn]+=internal_data_pt(0)->value(0)-
450  (*Target_rotation_angle_pt);
451 #ifdef PARANOID
452  }
453  else
454  {
455  std::stringstream junk;
456  junk << "Target_rotation_angle_pt=0\n"
457  << "Set with impose_rotation_angle(...)\n";
458  throw OomphLibError(
459  junk.str(),
462  }
463 #endif
464  }
465  }
466  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ impose_horizontal_force() [1/2]

void CircularPenetratorElement::impose_horizontal_force ( double target_horizontal_force_pt)
inline

Impose horizontal force (rather than rotation about contact node). Target force specified via pointer.

550  {
551  Target_horizontal_force_pt=target_horizontal_force_pt;
553  }

◆ impose_horizontal_force() [2/2]

void CircularPenetratorElement::impose_horizontal_force ( double target_horizontal_force_pt)
inline

Impose horizontal force (rather than rotation about contact node). Target force specified via pointer.

374  {
375  Target_horizontal_force_pt=target_horizontal_force_pt;
377  }

◆ impose_rotation_angle() [1/2]

void CircularPenetratorElement::impose_rotation_angle ( double target_rotation_angle_pt)
inline

Impose rotation about contact node (rather than horizontal force) Target angle specified via pointer.

558  {
560  Target_rotation_angle_pt=target_rotation_angle_pt;
561  }

◆ impose_rotation_angle() [2/2]

void CircularPenetratorElement::impose_rotation_angle ( double target_rotation_angle_pt)
inline

Impose rotation about contact node (rather than horizontal force) Target angle specified via pointer.

382  {
384  Target_rotation_angle_pt=target_rotation_angle_pt;
385  }

◆ impose_weight() [1/2]

void CircularPenetratorElement::impose_weight ( double target_weight_pt)
inline

Impose weight (rather than imposed displacement). Target weight specified via pointer.

525  {
526  Target_weight_pt=target_weight_pt;
527  Target_yc_pt=0;
528  }

◆ impose_weight() [2/2]

void CircularPenetratorElement::impose_weight ( double target_weight_pt)
inline

Impose weight (rather than imposed displacement). Target weight specified via pointer.

349  {
350  Target_weight_pt=target_weight_pt;
351  Target_yc_pt=0;
352  }

◆ impose_yc() [1/2]

void CircularPenetratorElement::impose_yc ( double target_yc_pt)
inline

Impose vertical position of control node (rather than weight). Target vertical position of control node specified via pointer.

533  {
535  Target_yc_pt=target_yc_pt;
536  }

◆ impose_yc() [2/2]

void CircularPenetratorElement::impose_yc ( double target_yc_pt)
inline

Impose vertical position of control node (rather than weight). Target vertical position of control node specified via pointer.

357  {
359  Target_yc_pt=target_yc_pt;
360  }

◆ output() [1/2]

void CircularPenetratorElement::output ( std::ostream &  outfile,
const unsigned nplot 
) const
inlinevirtual

Output coordinates of penetrator at nplot plot points.

Implements oomph::Penetrator.

721  {
722  for (unsigned j=0;j<nplot;j++)
723  {
724  double phi=2.0*MathematicalConstants::Pi*double(j)/double(nplot-1);
725  outfile << centre(0)+(*Radius_pt)*cos(phi) << " "
726  << centre(1)+(*Radius_pt)*sin(phi)
727  << std::endl;
728  }
729  }
double Pi
Definition: two_d_biharmonic.cc:235
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References cos(), j, BiharmonicTestFunctions2::Pi, and sin().

◆ output() [2/2]

void CircularPenetratorElement::output ( std::ostream &  outfile,
const unsigned nplot 
) const
inlinevirtual

Output coordinates of penetrator at nplot plot points.

Implements oomph::Penetrator.

539  {
540  for (unsigned j=0;j<nplot;j++)
541  {
542  double phi=2.0*MathematicalConstants::Pi*double(j)/double(nplot-1);
543  outfile << centre(0)+(*Radius_pt)*cos(phi) << " "
544  << centre(1)+(*Radius_pt)*sin(phi)
545  << std::endl;
546  }
547  }

References cos(), j, BiharmonicTestFunctions2::Pi, and sin().

◆ penetration() [1/2]

void CircularPenetratorElement::penetration ( const Vector< double > &  x,
const Vector< double > &  n,
double d,
bool intersection 
) const
inlinevirtual

Get penetration for given point x.

Implements oomph::Penetrator.

687  {
688 
689  // Vector from potential contact point to centre of penetrator
690  Vector<double> l(2);
691  l[0]=centre(0)-x[0];
692  l[1]=centre(1)-x[1];
693 
694  // Distance from potential contact point to centre of penetrator
695  double ll=sqrt(l[0]*l[0]+l[1]*l[1]);
696 
697  // Projection of vector from potential contact point to centre of penetrator
698  // onto outer unit normal on potential contact point
699  double project=n[0]*l[0]+n[1]*l[1];
700  double project_squared=project*project;
701 
702  // Final term in square root
703  double b_squared=ll*ll-(*Radius_pt)*(*Radius_pt);
704 
705  // Is square root negative? In this case we have no intersection
706  if (project_squared<b_squared)
707  {
708  d = -DBL_MAX;
709  intersection = false;
710  }
711  else
712  {
713  double sqr=sqrt(project_squared-b_squared);
714  d = -std::min(project-sqr,project+sqr);
715  intersection = true;
716  }
717  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define min(a, b)
Definition: datatypes.h:22
list x
Definition: plotDoE.py:28

References min, n, sqrt(), and plotDoE::x.

◆ penetration() [2/2]

void CircularPenetratorElement::penetration ( const Vector< double > &  x,
const Vector< double > &  n,
double d,
bool intersection 
) const
inlinevirtual

Get penetration for given point x.

Implements oomph::Penetrator.

504  {
505  // Vector from potential contact point to centre of penetrator
506  Vector<double> l(2);
507  l[0]=centre(0)-x[0];
508  l[1]=centre(1)-x[1];
509 
510  // Distance from potential contact point to centre of penetrator
511  double ll=sqrt(l[0]*l[0]+l[1]*l[1]);
512 
513  // Projection of vector from potential contact point to centre of penetrator
514  // onto outer unit normal on potential contact point
515  double project=n[0]*l[0]+n[1]*l[1];
516  double project_squared=project*project;
517 
518  // Final term in square root
519  double b_squared=ll*ll-(*Radius_pt)*(*Radius_pt);
520 
521  // Is square root negative? In this case we have no intersection
522  // and we return penetration as -DBL_MAX
523  if (project_squared<b_squared)
524  {
525  d= -DBL_MAX;
526  intersection = false;
527  }
528  else
529  {
530  double sqr=sqrt(project_squared-b_squared);
531  d= -std::min(project-sqr,project+sqr);
532  intersection = true;
533  }
534  }

References min, n, sqrt(), and plotDoE::x.

◆ radius() [1/2]

double CircularPenetratorElement::radius ( ) const
inline

Radius of penetrator.

751 {return *Radius_pt;}

◆ radius() [2/2]

double CircularPenetratorElement::radius ( ) const
inline

Radius of penetrator.

569 {return *Radius_pt;}

◆ resulting_force() [1/2]

Vector<double> CircularPenetratorElement::resulting_force ( ) const
inline

Resulting force from all associated ContactElements.

733  {
734  Vector<double> contact_force(2,0.0);
735  Vector<double> el_contact_force(2);
736  unsigned nel=Contact_element_mesh_pt->nelement();
737  for (unsigned e=0;e<nel;e++)
738  {
739  dynamic_cast<TemplateFreeContactElementBase*>(
741  resulting_contact_force(el_contact_force);
742  for (unsigned i=0;i<2;i++)
743  {
744  contact_force[i]+=el_contact_force[i];
745  }
746  }
747  return contact_force;
748  }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
Template-free base class for contact elements.
Definition: contact_elements.h:518

References e(), and i.

◆ resulting_force() [2/2]

Vector<double> CircularPenetratorElement::resulting_force ( ) const
inline

Resulting force from all associated ContactElements.

551  {
552  Vector<double> contact_force(2,0.0);
553  Vector<double> el_contact_force(2);
554  unsigned nel=Contact_element_mesh_pt->nelement();
555  for (unsigned e=0;e<nel;e++)
556  {
557  dynamic_cast<TemplateFreeContactElementBase*>(
559  resulting_contact_force(el_contact_force);
560  for (unsigned i=0;i<2;i++)
561  {
562  contact_force[i]+=el_contact_force[i];
563  }
564  }
565  return contact_force;
566  }

References e(), and i.

◆ rotation_angle_is_imposed() [1/2]

bool CircularPenetratorElement::rotation_angle_is_imposed ( )
inline

Is angle of rotation about control node imposed? If false then horizontal force is imposed.

542  {
543  return (Target_rotation_angle_pt!=0);
544  }

Referenced by ContactProblem< ELEMENT >::complete_problem_setup().

◆ rotation_angle_is_imposed() [2/2]

bool CircularPenetratorElement::rotation_angle_is_imposed ( )
inline

Is angle of rotation about control node imposed? If false then horizontal force is imposed.

366  {
367  return (Target_rotation_angle_pt!=0);
368  }

◆ set_angle() [1/2]

void CircularPenetratorElement::set_angle ( const double angle)
inline

Set angle of rotation around contact point.

362  {
364  }

References Jeffery_Solution::angle().

◆ set_angle() [2/2]

void CircularPenetratorElement::set_angle ( const double angle)
inline

Set angle of rotation around contact point.

193  {
195  }

References Jeffery_Solution::angle().

◆ set_contact_element_mesh_pt() [1/2]

void CircularPenetratorElement::set_contact_element_mesh_pt ( Mesh contact_element_mesh_pt)
inline

Set pointer to mesh of contact elements and setup external Data, i.e. Data that affects the residuals in this element. Also set the node pointed to by Control_node_pt as external Data for the elements in the contact mesh (unless they contain this node already).

380  {
383 
384  // Store Data associated with control node: It contains the traded
385  // Lagrange multiplier (contact pressure)
388 
389  // Store its position data
391 
392  // Store it as Data (which includes the displacement)
394 
395  // Loop over all the elements in the contact mesh
396  // If they don't contain the contact node already, its position
397  // is external data because it affects the penetration.
398  unsigned nel=Contact_element_mesh_pt->nelement();
399  for (unsigned e=0;e<nel;e++)
400  {
401  bool el_contains_control_node=false;
403  unsigned nnod=el_pt->nnode();
404  for (unsigned j=0;j<nnod;j++)
405  {
406  SolidNode* nod_pt=dynamic_cast<SolidNode*>(el_pt->node_pt(j));
407  if (nod_pt==Control_node_pt)
408  {
409  el_contains_control_node=true;
410  }
411  }
412 
413  // Position of control node affects position of penetrator and
414  // therefore is external data for all contact elements (apart from
415  // any that contain the control node as one of their own)
416  if (!el_contains_control_node)
417  {
418  // position
420  // displacement relative to position
422  }
423 
424  // Rotation angle angle affects position of penetrator and therefore
425  // affects penetration at all contact elements
427  }
428  }
Mesh * contact_element_mesh_pt() const
Definition: linear_solid_contact_with_gravity.cc:369
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void flush_external_data()
Flush all external data.
Definition: elements.cc:387
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:307
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
Definition: nodes.h:1686
Data *const & variable_position_pt() const
Pointer to variable_position data (const version)
Definition: nodes.h:1765

References oomph::GeneralisedElement::add_external_data(), e(), j, oomph::FiniteElement::nnode(), and oomph::FiniteElement::node_pt().

Referenced by ContactProblem< ELEMENT >::complete_problem_setup().

◆ set_contact_element_mesh_pt() [2/2]

void CircularPenetratorElement::set_contact_element_mesh_pt ( Mesh contact_element_mesh_pt)
inline

Set pointer to mesh of contact elements and setup external Data, i.e. Data that affects the residuals in this element. Also set the node pointed to by Control_node_pt as external Data for the elements in the contact mesh (unless they contain this node already).

211  {
214 
215  // Store Data associated with control node: It contains the traded
216  // Lagrange multiplier (contact pressure)
219 
220  // Store its position data
222 
223  // Loop over all the elements in the contact mesh
224  // If they don't contain the contact node already, its position
225  // is external data because it affects the penetration.
226  unsigned nel=Contact_element_mesh_pt->nelement();
227  for (unsigned e=0;e<nel;e++)
228  {
229  bool el_contains_control_node=false;
231  unsigned nnod=el_pt->nnode();
232  for (unsigned j=0;j<nnod;j++)
233  {
234  SolidNode* nod_pt=dynamic_cast<SolidNode*>(el_pt->node_pt(j));
235  if (nod_pt==Control_node_pt)
236  {
237  el_contains_control_node=true;
238  }
239  }
240 
241  // Position of control node affects position of penetrator and
242  // therefore is external data for all contact elements (apart from
243  // any that contain the control node as one of their own)
244  if (!el_contains_control_node)
245  {
247  }
248 
249  // Rotation angle angle affects position of penetrator and therefore
250  // affects penetration at all contact elements
252  }
253  }

References oomph::GeneralisedElement::add_external_data(), e(), j, oomph::FiniteElement::nnode(), and oomph::FiniteElement::node_pt().

◆ set_equilibrium_target_forces() [1/2]

void CircularPenetratorElement::set_equilibrium_target_forces ( )
inline

Set target horizontal and vertical force to be in equilibrium.

435  {
436 #ifdef PARANOID
438  {
439  std::stringstream junk;
440  junk << "Target_horizontal_force_pt==0\n"
441  << "Please set it by call to impose_weight(...)\n";
442  throw OomphLibError(
443  junk.str(),
446  }
447  if (Target_weight_pt==0)
448  {
449  std::stringstream junk;
450  junk << "Target_weight_pt==0.\n"
451  << "Please set it by call to impose_horizontal_force(...)\n";
452  throw OomphLibError(
453  junk.str(),
456  }
457 #endif
459  (*Target_horizontal_force_pt)=-force[0];
460  (*Target_weight_pt)=-force[1];
461  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ set_equilibrium_target_forces() [2/2]

void CircularPenetratorElement::set_equilibrium_target_forces ( )
inline

Set target horizontal and vertical force to be in equilibrium.

260  {
261 #ifdef PARANOID
263  {
264  std::stringstream junk;
265  junk << "Target_horizontal_force_pt==0\n"
266  << "Please set it by call to impose_weight(...)\n";
267  throw OomphLibError(
268  junk.str(),
271  }
272  if (Target_weight_pt==0)
273  {
274  std::stringstream junk;
275  junk << "Target_weight_pt==0.\n"
276  << "Please set it by call to impose_horizontal_force(...)\n";
277  throw OomphLibError(
278  junk.str(),
281  }
282 #endif
284  (*Target_horizontal_force_pt)=-force[0];
285  (*Target_weight_pt)=-force[1];
286  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ target_horizontal_force() [1/2]

double CircularPenetratorElement::target_horizontal_force ( )
inline

Target horizontal force (returns zero if not imposed)

479  {
481  {
482  return 0.0;
483  }
484  else
485  {
487  }
488  }

◆ target_horizontal_force() [2/2]

double CircularPenetratorElement::target_horizontal_force ( )
inline

Target horizontal force (returns zero if not imposed)

303  {
305  {
306  return 0.0;
307  }
308  else
309  {
311  }
312  }

◆ target_rotation_angle() [1/2]

double CircularPenetratorElement::target_rotation_angle ( )
inline

Target rotation angle about contact point (returns zero if not imposed)

505  {
507  {
508  return 0.0;
509  }
510  else
511  {
512  return *Target_rotation_angle_pt;
513  }
514  }

◆ target_rotation_angle() [2/2]

double CircularPenetratorElement::target_rotation_angle ( )
inline

Target rotation angle about contact point (returns zero if not imposed)

329  {
331  {
332  return 0.0;
333  }
334  else
335  {
336  return *Target_rotation_angle_pt;
337  }
338  }

◆ target_weight() [1/2]

double CircularPenetratorElement::target_weight ( )
inline

Target weight (returns zero if not imposed)

466  {
467  if (Target_weight_pt==0)
468  {
469  return 0.0;
470  }
471  else
472  {
473  return *Target_weight_pt;
474  }
475  }

◆ target_weight() [2/2]

double CircularPenetratorElement::target_weight ( )
inline

Target weight (returns zero if not imposed)

290  {
291  if (Target_weight_pt==0)
292  {
293  return 0.0;
294  }
295  else
296  {
297  return *Target_weight_pt;
298  }
299  }

◆ target_yc() [1/2]

double CircularPenetratorElement::target_yc ( )
inline

Target vertical position of control point (returns zero if not imposed)

492  {
493  if (Target_yc_pt==0)
494  {
495  return 0.0;
496  }
497  else
498  {
499  return *Target_yc_pt;
500  }
501  }

◆ target_yc() [2/2]

double CircularPenetratorElement::target_yc ( )
inline

Target vertical position of control point (returns zero if not imposed)

316  {
317  if (Target_yc_pt==0)
318  {
319  return 0.0;
320  }
321  else
322  {
323  return *Target_yc_pt;
324  }
325  }

◆ yc_is_imposed() [1/2]

bool CircularPenetratorElement::yc_is_imposed ( )
inline

Is vertical positon of control node imposed? If false then weight imposed.

518  {
519  return (Target_yc_pt!=0);
520  }

Referenced by ContactProblem< ELEMENT >::complete_problem_setup().

◆ yc_is_imposed() [2/2]

bool CircularPenetratorElement::yc_is_imposed ( )
inline

Is vertical positon of control node imposed? If false then weight imposed.

342  {
343  return (Target_yc_pt!=0);
344  }

Member Data Documentation

◆ Contact_element_mesh_pt

Mesh * CircularPenetratorElement::Contact_element_mesh_pt
private

Mesh of contact elements that contribute to weight/horizontal force.

◆ Control_node_pt

SolidNode * CircularPenetratorElement::Control_node_pt
private

Control node.

◆ External_data_index_of_traded_contact_pressure

unsigned CircularPenetratorElement::External_data_index_of_traded_contact_pressure
private

Index of external data that contains the the contact pressure in its Index_of_contact_pressure-th value

◆ Index_of_contact_pressure

unsigned CircularPenetratorElement::Index_of_contact_pressure
private

Index at which contact pressure (Lagr mult) is stored in nodal data associated with control node

◆ Index_of_horizontal_displacement

unsigned CircularPenetratorElement::Index_of_horizontal_displacement
private

Where is the horizontal displacement (linear_elasticity) stored?

◆ Index_of_vertical_displacement

unsigned CircularPenetratorElement::Index_of_vertical_displacement
private

Where is the vertical displacement (linear_elasticity) stored?

◆ Radius_pt

double * CircularPenetratorElement::Radius_pt
private

Pointer to radius of penetrator.

◆ Target_horizontal_force_pt

double * CircularPenetratorElement::Target_horizontal_force_pt
private

Pointer to target horizontal force (null if rotation angle angle about control node is imposed)

◆ Target_rotation_angle_pt

double * CircularPenetratorElement::Target_rotation_angle_pt
private

Pointer to target rotation angle about control node (null if horizontal force is imposed)

◆ Target_weight_pt

double * CircularPenetratorElement::Target_weight_pt
private

Pointer to target weight (null if vertical displacement of control node is imposed)

◆ Target_yc_pt

double * CircularPenetratorElement::Target_yc_pt
private

Pointer to target vertical displacement of control node (null if weight is imposed)


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