oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID > Class Template Reference

#include <pseudosolid_node_update_elements.h>

+ Inheritance diagram for oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >:

Public Member Functions

 PseudoSolidNodeUpdateElement ()
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void compute_norm (double &el_norm)
 
unsigned required_nvalue (const unsigned &n) const
 The required number of values is the sum of the two. More...
 
int solid_p_nodal_index () const
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
void evaluate_shape_derivs_by_direct_fd ()
 Evaluate shape derivatives by direct finite differencing. More...
 
void evaluate_shape_derivs_by_chain_rule ()
 Evaluate shape derivatives by chain rule. More...
 
void fill_in_shape_derivatives (DenseMatrix< double > &jacobian)
 
void fill_in_shape_derivatives_by_fd (DenseMatrix< double > &jacobian)
 
void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
void output (std::ostream &outfile)
 Overload the output function: Call that of the basic element. More...
 
void output (std::ostream &outfile, const unsigned &n_p)
 
void output (FILE *file_pt)
 Overload the output function: Call that of the basic element. More...
 
void output (FILE *file_pt, const unsigned &n_p)
 Output function is just the same as the basic equations. More...
 
unsigned num_Z2_flux_terms ()
 
void compute_exact_Z2_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm)
 
void get_Z2_flux (const Vector< double > &s, Vector< double > &flux)
 
unsigned nvertex_node () const
 Number of vertex nodes in the element. More...
 
Nodevertex_node_pt (const unsigned &j) const
 Pointer to the j-th vertex node in the element. More...
 
unsigned nrecovery_order ()
 
unsigned ndof_types () const
 
unsigned nbasic_dof_types () const
 
unsigned nsolid_dof_types () const
 
void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 

Private Attributes

bool Shape_derivs_by_direct_fd
 Boolean flag to indicate shape derivative method. More...
 

Detailed Description

template<class BASIC, class SOLID>
class oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >

A templated class that permits combination two different element types, for the solution of problems in deforming domains. The first template paremter BASIC is the standard element and the second SOLID solves the equations that are used to control the mesh deformation.

Constructor & Destructor Documentation

◆ PseudoSolidNodeUpdateElement()

template<class BASIC , class SOLID >
oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::PseudoSolidNodeUpdateElement ( )
inline

Constructor, call the BASIC and SOLID elements' constructors and set the "density" parameter for solid element to zero

67  {
68  SOLID::lambda_sq_pt() = &PseudoSolidHelper::Zero;
69  }
bool Shape_derivs_by_direct_fd
Boolean flag to indicate shape derivative method.
Definition: pseudosolid_node_update_elements.h:60
double Zero
Definition: pseudosolid_node_update_elements.cc:35

References oomph::PseudoSolidHelper::Zero.

Member Function Documentation

◆ compute_exact_Z2_error()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::compute_exact_Z2_error ( std::ostream &  outfile,
FiniteElement::SteadyExactSolutionFctPt  exact_flux_pt,
double error,
double norm 
)
inline

Plot the error when compared against a given exact flux. Also calculates the norm of the error and that of the exact flux. Use version in BASIC element

468  {
469  BASIC::compute_exact_Z2_error(outfile, exact_flux_pt, error, norm);
470  }
int error
Definition: calibrate.py:297

References calibrate::error.

◆ compute_norm()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::compute_norm ( double el_norm)
inline

Compute norm of solution: use the version in the BASIC class if there's any ambiguity

89  {
90  BASIC::compute_norm(el_norm);
91  }

◆ describe_local_dofs()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::describe_local_dofs ( std::ostream &  out,
const std::string &  current_string 
) const
inline

Function to describe the local dofs of the element. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...)

81  {
82  BASIC::describe_local_dofs(out, current_string);
83  SOLID::describe_local_dofs(out, current_string);
84  }
std::ofstream out("Result.txt")

References out().

◆ evaluate_shape_derivs_by_chain_rule()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::evaluate_shape_derivs_by_chain_rule ( )
inline

Evaluate shape derivatives by chain rule.

168  {
170  }

References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::Shape_derivs_by_direct_fd.

◆ evaluate_shape_derivs_by_direct_fd()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::evaluate_shape_derivs_by_direct_fd ( )
inline

Evaluate shape derivatives by direct finite differencing.

162  {
164  }

References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::Shape_derivs_by_direct_fd.

◆ fill_in_contribution_to_jacobian()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inline

Final override for jacobian function: Contributions are included from both the underlying element types

131  {
132  // Call the basic equations first
133  BASIC::fill_in_contribution_to_jacobian(residuals, jacobian);
134  // Call the solid equations
135  SOLID::fill_in_contribution_to_jacobian(residuals, jacobian);
136 
137  // Now fill in the off-diagonal entries (the shape derivatives),
138  fill_in_shape_derivatives(jacobian);
139  }
void fill_in_shape_derivatives(DenseMatrix< double > &jacobian)
Definition: pseudosolid_node_update_elements.h:175
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: gen_axisym_advection_diffusion_elements.h:536

References oomph::fill_in_contribution_to_jacobian(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().

◆ fill_in_contribution_to_jacobian_and_mass_matrix()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_jacobian_and_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix 
)
inline

Final override for mass matrix function: contributions are included from both the underlying element types

147  {
148  // Call the basic equations first
150  residuals, jacobian, mass_matrix);
151  // Call the solid equations
153  residuals, jacobian, mass_matrix);
154 
155  // Now fill in the off-diagonal entries (the shape derivatives),
156  fill_in_shape_derivatives(jacobian);
157  }
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: gen_axisym_advection_diffusion_elements.h:547

References oomph::fill_in_contribution_to_jacobian_and_mass_matrix(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().

◆ fill_in_contribution_to_residuals()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inline

Final override for the residuals function. Contributions are added from both underlying element types

120  {
121  // Call the basic equations first
123  // Add the solid equations contribution
125  }
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: gen_axisym_advection_diffusion_elements.h:522

References oomph::fill_in_contribution_to_residuals().

◆ fill_in_shape_derivatives()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives ( DenseMatrix< double > &  jacobian)
inline

Fill in the shape derivatives of the BASIC equations w.r.t. the solid position dofs

176  {
177  // Default is to use finite differences
179  {
180  this->fill_in_shape_derivatives_by_fd(jacobian);
181  }
182  // Otherwise need to do a bit more work
183  else
184  {
185  // Calculate storage requirements
186  const unsigned n_dof = this->ndof();
187  const unsigned n_node = this->nnode();
188  const unsigned nodal_dim = this->nodal_dimension();
189 
190  // If there are no nodes or dofs return
191  if ((n_dof == 0) || (n_node == 0))
192  {
193  return;
194  }
195 
196  // Generalised dofs have NOT been considered, shout
197  if (this->nnodal_position_type() != 1)
198  {
199  throw OomphLibError("Shape derivatives do not (yet) allow for "
200  "generalised position dofs\n",
203  }
204 
205  // Storage for derivatives of residuals w.r.t. nodal coordinates
206  RankThreeTensor<double> dresidual_dnodal_coordinates(
207  n_dof, nodal_dim, n_node, 0.0);
208 
209  // Get the analytic derivatives for the BASIC equations
210  BASIC::get_dresidual_dnodal_coordinates(dresidual_dnodal_coordinates);
211 
212  // Now add the appropriate contributions to the Jacobian
213  int local_unknown = 0;
214 
215  // Loop over dofs
216  //(this will include the solid dofs,
217  // but all those contributions should be zero)
218  for (unsigned l = 0; l < n_dof; l++)
219  {
220  // Loop over the nodes
221  for (unsigned n = 0; n < n_node; n++)
222  {
223  // Loop over the position_types (only one)
224  unsigned k = 0;
225  // Loop over the coordinates
226  for (unsigned i = 0; i < nodal_dim; i++)
227  {
228  // Get the equation of the local unknown
229  local_unknown = this->position_local_eqn(n, k, i);
230 
231  // If not pinned, add the contribution to the Jacobian
232  if (local_unknown >= 0)
233  {
234  jacobian(l, local_unknown) +=
235  dresidual_dnodal_coordinates(l, i, n);
236  }
237  }
238  }
239  }
240  }
241  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
void fill_in_shape_derivatives_by_fd(DenseMatrix< double > &jacobian)
Definition: pseudosolid_node_update_elements.h:246
char char char int int * k
Definition: level2_impl.h:374
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd(), i, k, n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::Shape_derivs_by_direct_fd.

Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_jacobian(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_jacobian_and_mass_matrix().

◆ fill_in_shape_derivatives_by_fd()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd ( DenseMatrix< double > &  jacobian)
inline

Fill in the derivatives of the BASIC equations w.r.t. the solid position dofs

247  {
248  // Flag to indicate if we use first or second order FD
249  // bool use_first_order_fd=false;
250 
251  // Find the number of nodes
252  const unsigned n_node = this->nnode();
253 
254  // If there aren't any nodes, then return straight away
255  if (n_node == 0)
256  {
257  return;
258  }
259 
260  // Call the update function to ensure that the element is in
261  // a consistent state before finite differencing starts
262  this->update_before_solid_position_fd();
263 
264  // Get the number of position dofs and dimensions at the node
265  const unsigned n_position_type = this->nnodal_position_type();
266  const unsigned nodal_dim = this->nodal_dimension();
267 
268  // Find the number of dofs in the element
269  const unsigned n_dof = this->ndof();
270 
271  // Create residual newres vectors
272  Vector<double> residuals(n_dof);
273  Vector<double> newres(n_dof);
274  // Vector<double> newres_minus(n_dof);
275 
276  // Calculate the residuals (for the BASIC) equations
277  // Need to do this using fill_in because get_residuals will
278  // compute all residuals for the problem, which is
279  // a little ineffecient
280  for (unsigned m = 0; m < n_dof; m++)
281  {
282  residuals[m] = 0.0;
283  }
285 
286  // Need to determine which degrees of freedom are solid degrees of
287  // freedom
288  // A vector of booleans that will be true if the dof is associated
289  // with the solid equations
290  std::vector<bool> dof_is_solid(n_dof, false);
291 
292  // Now set all solid positional dofs in the vector
293  for (unsigned n = 0; n < n_node; n++)
294  {
295  for (unsigned k = 0; k < n_position_type; k++)
296  {
297  for (unsigned i = 0; i < nodal_dim; i++)
298  {
299  int local_dof = this->position_local_eqn(n, k, i);
300  if (local_dof >= 0)
301  {
302  dof_is_solid[local_dof] = true;
303  }
304  }
305  }
306  }
307 
308  // Add the solid pressures (in solid elements without
309  // solid pressure the number will be zero).
310  unsigned n_solid_pres = this->npres_solid();
311  for (unsigned l = 0; l < n_solid_pres; l++)
312  {
313  int local_dof = this->solid_p_local_eqn(l);
314  if (local_dof >= 0)
315  {
316  dof_is_solid[local_dof] = true;
317  }
318  }
319 
320 
321  // Integer storage for local unknown
322  int local_unknown = 0;
323 
324  // Use default value defined in GeneralisedElement
325  const double fd_step = this->Default_fd_jacobian_step;
326 
327  // Loop over the nodes
328  for (unsigned n = 0; n < n_node; n++)
329  {
330  // Loop over position dofs
331  for (unsigned k = 0; k < n_position_type; k++)
332  {
333  // Loop over dimension
334  for (unsigned i = 0; i < nodal_dim; i++)
335  {
336  // If the variable is free
337  local_unknown = this->position_local_eqn(n, k, i);
338  if (local_unknown >= 0)
339  {
340  // Store a pointer to the (generalised) Eulerian nodal position
341  double* const value_pt = &(this->node_pt(n)->x_gen(k, i));
342 
343  // Save the old value of the (generalised) Eulerian nodal position
344  const double old_var = *value_pt;
345 
346  // Increment the (generalised) Eulerian nodal position
347  *value_pt += fd_step;
348 
349  // Perform any auxialiary node updates
350  this->node_pt(n)->perform_auxiliary_node_update_fct();
351 
352  // Calculate the new residuals
353  // Need to do this using fill_in because get_residuals will
354  // compute all residuals for the problem, which is
355  // a little ineffecient
356  for (unsigned m = 0; m < n_dof; m++)
357  {
358  newres[m] = 0.0;
359  }
361 
362  // if (use_first_order_fd)
363  {
364  // Do forward finite differences
365  for (unsigned m = 0; m < n_dof; m++)
366  {
367  // Stick the entry into the Jacobian matrix
368  // but only if it's not a solid dof
369  if (dof_is_solid[m] == false)
370  {
371  jacobian(m, local_unknown) =
372  (newres[m] - residuals[m]) / fd_step;
373  }
374  }
375  }
376  // else
377  // {
378  // //Take backwards step for the (generalised)
379  // Eulerian nodal
380  // // position
381  // node_pt(n)->x_gen(k,i) = old_var-fd_step;
382 
383  // //Calculate the new residuals at backward position
384  // //BASIC::get_residuals(newres_minus);
385 
386  // //Do central finite differences
387  // for(unsigned m=0;m<n_dof;m++)
388  // {
389  // //Stick the entry into the Jacobian matrix
390  // jacobian(m,local_unknown) =
391  // (newres[m] - newres_minus[m])/(2.0*fd_step);
392  // }
393  // }
394 
395  // Reset the (generalised) Eulerian nodal position
396  *value_pt = old_var;
397 
398  // Perform any auxialiary node updates
399  this->node_pt(n)->perform_auxiliary_node_update_fct();
400  }
401  }
402  }
403  }
404 
405  // End of finite difference loop
406  // Final reset of any dependent data
407  this->reset_after_solid_position_fd();
408  }
int * m
Definition: level2_cplx_impl.h:294

References oomph::fill_in_contribution_to_residuals(), i, k, m, and n.

Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().

◆ get_dof_numbers_for_unknowns()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::get_dof_numbers_for_unknowns ( std::list< std::pair< unsigned long, unsigned >> &  dof_lookup_list) const
inline

Create a list of pairs for all unknowns in this element, so that the first entry in each pair contains the global equation number of the unknown, while the second one contains the number of the "DOF type" that this unknown is associated with. This method combines the get_dof_numbers_for_unknowns(...) method for the BASIC and SOLID elements. The basic elements retain their DOF type numbering and the SOLID elements DOF type numbers are incremented by nbasic_dof_types().

530  {
531  // get the solid list
532  std::list<std::pair<unsigned long, unsigned>> solid_list;
533  SOLID::get_dof_numbers_for_unknowns(solid_list);
534 
535  // get the basic list
536  BASIC::get_dof_numbers_for_unknowns(dof_lookup_list);
537 
538  // get the number of basic dof types
539  unsigned nbasic_dof_types = BASIC::ndof_types();
540 
541  // add the solid lookup list to the basic lookup list
542  // incrementing the solid dof numbers by nbasic_dof_types
543  typedef std::list<std::pair<unsigned long, unsigned>>::iterator IT;
544  for (IT it = solid_list.begin(); it != solid_list.end(); it++)
545  {
546  std::pair<unsigned long, unsigned> new_pair;
547  new_pair.first = it->first;
548  new_pair.second = it->second + nbasic_dof_types;
549  dof_lookup_list.push_front(new_pair);
550  }
551  }
unsigned nbasic_dof_types() const
Definition: pseudosolid_node_update_elements.h:508

References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::nbasic_dof_types().

◆ get_Z2_flux()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::get_Z2_flux ( const Vector< double > &  s,
Vector< double > &  flux 
)
inline

'Flux' vector for Z2 error estimation: Error estimation is based on error in BASIC element

475  {
477  }
RealScalar s
Definition: level1_cplx_impl.h:130
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery.
Definition: overloaded_element_body.h:745

References ProblemParameters::flux(), get_Z2_flux(), and s.

◆ identify_geometric_data()

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::identify_geometric_data ( std::set< Data * > &  geometric_data_pt)
inline

Specify Data that affects the geometry of the element by adding the position Data to the set that's passed in. (This functionality is required in FSI problems; set is used to avoid double counting).

416  {
417  // Loop over the node update data and add to the set
418  const unsigned n_node = this->nnode();
419  for (unsigned j = 0; j < n_node; j++)
420  {
421  geometric_data_pt.insert(
422  dynamic_cast<SolidNode*>(this->node_pt(j))->variable_position_pt());
423  }
424  }
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References j.

◆ nbasic_dof_types()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::nbasic_dof_types ( ) const
inline

return the number of DOF types associated with the BASIC elements in this combined element

509  {
510  return BASIC::ndof_types();
511  }

Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::get_dof_numbers_for_unknowns().

◆ ndof_types()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::ndof_types ( ) const
inline

The number of "DOF types" that degrees of freedom in this element are sub-divided into.

502  {
503  return BASIC::ndof_types() + SOLID::ndof_types();
504  }

◆ nrecovery_order()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::nrecovery_order ( )
inline

Order of recovery shape functions for Z2 error estimation: Done for BASIC element since it determines the refinement

494  {
495  return BASIC::nrecovery_order();
496  }

◆ nsolid_dof_types()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::nsolid_dof_types ( ) const
inline

return the number of DOF types associated with the SOLID elements in this combined element

516  {
517  return SOLID::ndof_types();
518  }

◆ num_Z2_flux_terms()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::num_Z2_flux_terms ( )
inline

Number of 'flux' terms for Z2 error estimation: Error estimation is based on error in BASIC element

455  {
456  return BASIC::num_Z2_flux_terms();
457  }

◆ nvertex_node()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::nvertex_node ( ) const
inline

Number of vertex nodes in the element.

481  {
482  return BASIC::nvertex_node();
483  }

◆ output() [1/4]

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::output ( FILE *  file_pt)
inline

Overload the output function: Call that of the basic element.

442  {
443  BASIC::output(file_pt);
444  }
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490

References output().

◆ output() [2/4]

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::output ( FILE *  file_pt,
const unsigned n_p 
)
inline

Output function is just the same as the basic equations.

448  {
449  BASIC::output(file_pt, n_p);
450  }

References output().

◆ output() [3/4]

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::output ( std::ostream &  outfile)
inline

Overload the output function: Call that of the basic element.

429  {
430  BASIC::output(outfile);
431  }

References output().

◆ output() [4/4]

template<class BASIC , class SOLID >
void oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::output ( std::ostream &  outfile,
const unsigned n_p 
)
inline

Output function: Plot at n_p plot points using the basic element's output function

436  {
437  BASIC::output(outfile, n_p);
438  }

References output().

◆ required_nvalue()

template<class BASIC , class SOLID >
unsigned oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::required_nvalue ( const unsigned n) const
inline

The required number of values is the sum of the two.

95  {
96  return BASIC::required_nvalue(n) + SOLID::required_nvalue(n);
97  }

References n.

◆ solid_p_nodal_index()

template<class BASIC , class SOLID >
int oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::solid_p_nodal_index ( ) const
inline

We assume that the solid stuff is stored at the end of the nodes, i.e. its index is the number of continuously interplated values in the BASIC equations.

103  {
104  // At the moment, we can't handle this case in generality so throw an
105  // error if the solid pressure is stored at the nodes
106  if (SOLID::solid_p_nodal_index() >= 0)
107  {
108  throw OomphLibError("Cannot handle (non-refineable) continuous solid "
109  "pressure interpolation",
112  }
113 
114  return SOLID::solid_p_nodal_index();
115  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ vertex_node_pt()

template<class BASIC , class SOLID >
Node* oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::vertex_node_pt ( const unsigned j) const
inline

Pointer to the j-th vertex node in the element.

487  {
488  return BASIC::vertex_node_pt(j);
489  }

References j.

Member Data Documentation

◆ Shape_derivs_by_direct_fd


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