oomph::SCoupledElement< ELEMENT > Class Template Reference

#include <SCoupledElement.h>

+ Inheritance diagram for oomph::SCoupledElement< ELEMENT >:

Public Member Functions

 SCoupledElement ()
 Constructor: Call constructor of underlying element. More...
 
 ~SCoupledElement ()
 Destructor (empty) More...
 
Node * construct_node (const unsigned &n)
 
Node * construct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
Node * construct_boundary_node (const unsigned &n)
 
Node * construct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Add the element's contribution to its residual vector (wrapper) More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Add the element's contribution to its residual vector and Jacobian matrix (wrapper) More...
 
void get_momentum_and_energy (double &mass, Vector< double > &lin_mo, Vector< double > &ang_mo, double &pot_en, double &kin_en)
 
void set_nodal_coupling_residual (const bool &isCoupled, Vector< Vector< double >> &residual)
 
double get_nodal_coupling_residual (const unsigned n, const unsigned d)
 

Private Member Functions

void add_external_coupling_forces_to_residuals (Vector< double > &residuals)
 Add the point source contribution to the residual vector. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 

Private Attributes

Vector< Vector< double > > nodal_coupling_residual
 Nodal coupling forces. More...
 

Detailed Description

template<class ELEMENT>
class oomph::SCoupledElement< ELEMENT >

Wrapper class for solid elements to be coupled with discrete solid particles on the surfaces

Constructor & Destructor Documentation

◆ SCoupledElement()

template<class ELEMENT >
oomph::SCoupledElement< ELEMENT >::SCoupledElement ( )
inline

Constructor: Call constructor of underlying element.

45  {};

◆ ~SCoupledElement()

template<class ELEMENT >
oomph::SCoupledElement< ELEMENT >::~SCoupledElement ( )
inline

Destructor (empty)

49  {};

Member Function Documentation

◆ add_external_coupling_forces_to_residuals()

template<class ELEMENT >
void oomph::SCoupledElement< ELEMENT >::add_external_coupling_forces_to_residuals ( Vector< double > &  residuals)
inlineprivate

Add the point source contribution to the residual vector.

300  {
301  // No further action if no coupling force
302  if (nodal_coupling_residual.size() == 0) return;
303 
304  // Find out how many nodes there are
305  const unsigned n_node = this->nnode();
306  const unsigned DIM = this->dim();
307 
308  //find out how many positional dofs there are
309  unsigned n_position_type = this->nnodal_position_type();
310 
311  // Set up memory for the shape/test functions
312  //Shape psi(n_node);
313 
314  // Loop over all nodes belonging to the element
315  int local_eqn = 0;
316  // Loop over the nodes of the element
317  for (unsigned l = 0; l < n_node; l++)
318  {
319  // Loop of types of dofs
320  for (unsigned k = 0; k < n_position_type; k++)
321  {
322  // Loop over the force components
323  for (unsigned i = 0; i < DIM; i++)
324  {
325  // Get the local equation
326  local_eqn = this->position_local_eqn(l, k, i);
327  // IF it's not a boundary condition and the interpolated force is nonzero
328  if (local_eqn >= 0)
329  {
330  // add the nodal coupling residual to the global residual vector
331  residuals[local_eqn] += nodal_coupling_residual[l][i];
332  }
333  }
334  }
335  }
336  }
Vector< Vector< double > > nodal_coupling_residual
Nodal coupling forces.
Definition: SCoupledElement.h:384
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References constants::i, and oomph::SCoupledElement< ELEMENT >::nodal_coupling_residual.

Referenced by oomph::SCoupledElement< ELEMENT >::fill_in_contribution_to_jacobian(), and oomph::SCoupledElement< ELEMENT >::fill_in_contribution_to_residuals().

◆ construct_boundary_node() [1/2]

template<class ELEMENT >
Node* oomph::SCoupledElement< ELEMENT >::construct_boundary_node ( const unsigned &  n)
inline

Construct the local node n and return a pointer to it. in the case when it is a boundary node; that is it MAY be located on a Mesh boundary Needs to be overridden because a SurfaceCoupledElement uses CoupledSolidNode's

94  {
95  // Construct a solid node and assign it to the local node pointer vector.
96  // The dimension and number of values are taken from internal element data
97  // The number of solid pressure dofs are also taken from internal data
98  // The number of timesteps to be stored comes from the problem!
99  this->node_pt(n) = new BoundaryNode<CoupledSolidNode>(
100  this->lagrangian_dimension(),
101  this->nnodal_lagrangian_type(),
102  this->nodal_dimension(),
103  this->nnodal_position_type(),
104  this->required_nvalue(n));
105  // Now return a pointer to the node, so that the mesh can find it
106  return this->node_pt(n);
107  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:32

References n.

◆ construct_boundary_node() [2/2]

template<class ELEMENT >
Node* oomph::SCoupledElement< ELEMENT >::construct_boundary_node ( const unsigned &  n,
TimeStepper *const &  time_stepper_pt 
)
inline

Construct the local node n and return a pointer to it, in the case when the node MAY be located on a boundary. Additionally, create storage for ‘history’ values as required by timestepper Needs to be overridden because a SurfaceCoupledElement uses CoupledSolidNode's

116  {
117  // Construct a solid node and assign it to the local node pointer vector
118  // The dimension and number of values are taken from internal element data
119  // The number of solid pressure dofs are also taken from internal data
120  this->node_pt(n) = new BoundaryNode<CoupledSolidNode>(
121  time_stepper_pt,
122  this->lagrangian_dimension(),
123  this->nnodal_lagrangian_type(),
124  this->nodal_dimension(),
125  this->nnodal_position_type(),
126  this->required_nvalue(n));
127  // Now return a pointer to the node, so that the mesh can find it
128  return this->node_pt(n);
129  }

References n.

◆ construct_node() [1/2]

template<class ELEMENT >
Node* oomph::SCoupledElement< ELEMENT >::construct_node ( const unsigned &  n)
inline

Construct the local node n and return a pointer to it. Needs to be overridden because a SurfaceCoupledElement uses CoupledSolidNode's

54  {
55  // Construct a solid node and assign it to the local node pointer vector.
56  // The dimension and number of values are taken from internal element data
57  // The number of solid pressure dofs are also taken from internal data
58  // The number of timesteps to be stored comes from the problem!
59  this->node_pt(n) = new CoupledSolidNode(
60  this->lagrangian_dimension(),
61  this->nnodal_lagrangian_type(),
62  this->nodal_dimension(),
63  this->nnodal_position_type(),
64  this->required_nvalue(n));
65  // Now return a pointer to the node, so that the mesh can find it
66  return this->node_pt(n);
67  }

References n.

◆ construct_node() [2/2]

template<class ELEMENT >
Node* oomph::SCoupledElement< ELEMENT >::construct_node ( const unsigned &  n,
TimeStepper *const &  time_stepper_pt 
)
inline

Construct the local node n and return a pointer to it. Additionally, create storage for ‘history’ values as required by timestepper Needs to be overridden because a SurfaceCoupledElement uses CoupledSolidNode's

74  {
75  // Construct a solid node and assign it to the local node pointer vector
76  // The dimension and number of values are taken from internal element data
77  // The number of solid pressure dofs are also taken from internal data
78  this->node_pt(n) = new CoupledSolidNode(
79  time_stepper_pt,
80  this->lagrangian_dimension(),
81  this->nnodal_lagrangian_type(),
82  this->nodal_dimension(),
83  this->nnodal_position_type(),
84  this->required_nvalue(n));
85  // Now return a pointer to the node, so that the mesh can find it
86  return this->node_pt(n);
87  }

References n.

◆ fill_in_contribution_to_jacobian()

template<class ELEMENT >
void oomph::SCoupledElement< ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inline

Add the element's contribution to its residual vector and Jacobian matrix (wrapper)

146  {
147  //Call the generic routine with the flag set to 1
148  ELEMENT::fill_in_generic_contribution_to_residuals_pvd(residuals, jacobian, 1);
149 
150  // Add point source contribution
152  }
void add_external_coupling_forces_to_residuals(Vector< double > &residuals)
Add the point source contribution to the residual vector.
Definition: SCoupledElement.h:299

References oomph::SCoupledElement< ELEMENT >::add_external_coupling_forces_to_residuals().

◆ fill_in_contribution_to_residuals()

template<class ELEMENT >
void oomph::SCoupledElement< ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inline

Add the element's contribution to its residual vector (wrapper)

134  {
135  // Call the generic residuals function with flag set to 0 using a dummy matrix argument
136  ELEMENT::fill_in_generic_contribution_to_residuals_pvd(
137  residuals, GeneralisedElement::Dummy_matrix, 0);
138 
139  // Add point source contribution
141  }

References oomph::SCoupledElement< ELEMENT >::add_external_coupling_forces_to_residuals().

◆ get_momentum_and_energy()

template<class ELEMENT >
void oomph::SCoupledElement< ELEMENT >::get_momentum_and_energy ( double mass,
Vector< double > &  lin_mo,
Vector< double > &  ang_mo,
double pot_en,
double kin_en 
)
inline
156  {
157  const unsigned DIM = this->dim();
158  // Initialise mass
159  mass = 0;
160  // Initialise momentum
161  lin_mo.initialise(0);
162  ang_mo.initialise(0);
163  // Initialise energy
164  pot_en = 0;
165  kin_en = 0;
166 
167  //Set the value of n_intpt
168  unsigned n_intpt = this->integral_pt()->nweight();
169 
170  //Set the Vector to hold local coordinates
171  Vector<double> s(DIM);
172 
173  //Find out how many nodes there are
174  const unsigned n_node = this->nnode();
175 
176  //Find out how many positional dofs there are
177  const unsigned n_position_type = this->nnodal_position_type();
178 
179  //Set up memory for the shape functions
180  Shape psi(n_node, n_position_type);
181  DShape dpsidxi(n_node, n_position_type, DIM);
182 
183  // Timescale ratio (non-dim density)
184  double lambda_sq = this->lambda_sq();
185 
186  //Loop over the integration points
187  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
188  {
189  //Assign values of s
190  for (unsigned i = 0; i < DIM; i++)
191  { s[i] = this->integral_pt()->knot(ipt, i); }
192 
193  //Get the integral weight
194  double w = this->integral_pt()->weight(ipt);
195 
196  //Call the derivatives of the shape functions and get Jacobian
197  double J = this->dshape_lagrangian_at_knot(ipt, psi, dpsidxi);
198 
199  //Storage for Lagrangian coordinates and velocity (initialised to zero)
200  Vector<double> interpolated_xi(DIM, 0.0);
201  Vector<double> veloc(DIM, 0.0);
202 
203  //Calculate lagrangian coordinates
204  for (unsigned l = 0; l < n_node; l++)
205  {
206  //Loop over positional dofs
207  for (unsigned k = 0; k < n_position_type; k++)
208  {
209  //Loop over displacement components (deformed position)
210  for (unsigned i = 0; i < DIM; i++)
211  {
212  //Calculate the Lagrangian coordinates
213  interpolated_xi[i] += this->lagrangian_position_gen(l, k, i) * psi(l, k);
214 
215  //Calculate the velocity components (if unsteady solve)
216  if (this->Unsteady)
217  {
218  veloc[i] += this->dnodal_position_gen_dt(l, k, i) * psi(l, k);
219  }
220  }
221  }
222  }
223 
224  //Get isotropic growth factor gamma
225  double gamma = 1.0;
226  this->get_isotropic_growth(ipt, s, interpolated_xi, gamma);
227 
228  //Premultiply the undeformed volume ratio (from the isotropic
229  // growth), the integral weights, the coupling weights, and the Jacobian
230  double W = gamma * w * J;
231 
232  DenseMatrix<double> sigma(DIM, DIM);
233  DenseMatrix<double> strain(DIM, DIM);
234 
235  //Now calculate the stress tensor from the constitutive law
236  this->get_stress(s, sigma);
237 
238  // Add pre-stress
239  for (unsigned i = 0; i < DIM; i++)
240  {
241  for (unsigned j = 0; j < DIM; j++)
242  {
243  sigma(i, j) += this->prestress(i, j, interpolated_xi);
244  }
245  }
246 
247  //get the strain
248  this->get_strain(s, strain);
249 
250  // Initialise
251  double local_pot_en = 0;
252  double veloc_sq = 0;
253 
254  // Compute integrals
255  for (unsigned i = 0; i < DIM; i++)
256  {
257  for (unsigned j = 0; j < DIM; j++)
258  {
259  local_pot_en += sigma(i, j) * strain(i, j);
260  }
261  veloc_sq += veloc[i] * veloc[i];
262  }
263 
264  // Mass
265  mass += lambda_sq * W;
266  // Linear momentum and angular momentum
267  Vector<double> cross_product(DIM, 0);
268  VectorHelpers::cross(interpolated_xi, veloc, cross_product);
269  for (unsigned i = 0; i < DIM; i++)
270  {
271  lin_mo[i] += lambda_sq * veloc[i] * W;
272  ang_mo[i] += lambda_sq * cross_product[i] * W;
273  }
274  // Potential energy
275  pot_en += 0.5 * local_pot_en * W;
276  // Kinetic energy
277  kin_en += lambda_sq * 0.5 * veloc_sq * W;
278  }
279  }
Mdouble gamma(Mdouble gamma_in)
This is the gamma function returns the true value for the half integer value.
Definition: ExtendedMath.cc:137

References mathsFunc::gamma(), and constants::i.

◆ get_nodal_coupling_residual()

template<class ELEMENT >
double oomph::SCoupledElement< ELEMENT >::get_nodal_coupling_residual ( const unsigned  n,
const unsigned  d 
)
inline
289  {
290  if (nodal_coupling_residual.size()>n && nodal_coupling_residual[n].size()>d)
291  return nodal_coupling_residual[n][d];
292  else
293  return 0;
294  }

References n, and oomph::SCoupledElement< ELEMENT >::nodal_coupling_residual.

◆ output()

template<class ELEMENT >
void oomph::SCoupledElement< ELEMENT >::output ( std::ostream &  outfile,
const unsigned &  n_plot 
)
inlineprivate
339  {
340  const unsigned DIM = this->dim();
341  Vector<double> x(DIM);
342  Vector<double> dxdt(DIM);
343  Vector<double> s(DIM);
344 
345  // Tecplot header info
346  outfile << this->tecplot_zone_string(n_plot);
347 
348  // Loop over plot points
349  unsigned num_plot_points = this->nplot_points(n_plot);
350  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
351  {
352  // Get local coordinates of plot point
353  this->get_s_plot(iplot, n_plot, s);
354 
355  // Get Eulerian coordinates
356  this->interpolated_x(s, x);
357  SolidFiniteElement* el_pt = dynamic_cast<SolidFiniteElement*>(this);
358  el_pt->interpolated_dxdt(s, 1, dxdt);
359 
360  // Dummy integration point
361  unsigned ipt = 0;
362 
363  // Output x,y,z
364  for (unsigned i = 0; i < DIM; i++)
365  {
366  outfile << x[i] << " ";
367  }
368 
369  // Output velocity dnodal_position_gen_dt
370  for (unsigned i = 0; i < DIM; i++)
371  {
372  outfile << dxdt[i] << " ";
373  }
374 
375  outfile << std::endl;
376  }
377 
378  // Write tecplot footer (e.g. FE connectivity lists)
379  this->write_tecplot_zone_footer(outfile, n_plot);
380  outfile << std::endl;
381  }

References constants::i.

◆ set_nodal_coupling_residual()

template<class ELEMENT >
void oomph::SCoupledElement< ELEMENT >::set_nodal_coupling_residual ( const bool isCoupled,
Vector< Vector< double >> &  residual 
)
inline
282  {
283  if (isCoupled)
284  { nodal_coupling_residual = residual; }
285  else
286  { nodal_coupling_residual.clear(); }
287  }

References oomph::SCoupledElement< ELEMENT >::nodal_coupling_residual.

Member Data Documentation

◆ nodal_coupling_residual


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