oomph::VolumeCoupledElement< ELEMENT > Class Template Reference

#include <VCoupledElement.h>

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

Public Member Functions

 VolumeCoupledElement ()
 Constructor: Call constructor of underlying element. More...
 
 ~VolumeCoupledElement ()
 Destructor (empty) More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals) override
 Add the element's contribution to its residual vector (wrapper) More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) override
 Add the element's contribution to its residual vector and Jacobian matrix (wrapper) More...
 
void set_nodal_coupling_residual (Vector< Vector< double > > &cResidual)
 
doubleget_nodal_coupling_residual (const unsigned &l, const unsigned &i)
 
void set_nodal_coupling_jacobian (Vector< Vector< double > > &cJacobian)
 
doubleget_nodal_coupling_jacobian (const unsigned &l, const unsigned &ll)
 
void setCouplingStiffness (Vector< Vector< double > > &cStiffness)
 
doublegetCouplingStiffness (const unsigned &m, const unsigned &l)
 
void get_momentum_and_energy (double &mass, Vector< double > &lin_mo, Vector< double > &ang_mo, double &pot_en, double &kin_en)
 

Private Member Functions

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

Private Attributes

Vector< Vector< double > > nodal_coupling_residual
 Nodal coupling forces. More...
 
Vector< Vector< double > > nodal_coupling_jacobian
 Nodal coupling jacobian. More...
 
Vector< Vector< double > > couplingStiffness
 Coupling stiffness to discrete particles (FIXME: should be moved into DPMVCoupledElement) More...
 

Detailed Description

template<class ELEMENT>
class oomph::VolumeCoupledElement< ELEMENT >

Wrapper class for solid elements to be coupled with discrete solid particles over a partially overlapped volume

Constructor & Destructor Documentation

◆ VolumeCoupledElement()

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

Constructor: Call constructor of underlying element.

56  {};

◆ ~VolumeCoupledElement()

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

Destructor (empty)

60  {};

Member Function Documentation

◆ add_internal_coupling_forces_to_residuals()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::add_internal_coupling_forces_to_residuals ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
const unsigned &  flag 
)
inlineprivate

Add the point source contribution to the residual vector.

624  {
625  // No further action if no coupling forces
626  if (nodal_coupling_residual.size() == 0) return;
627 
628  // Find out how many nodes there are
629  const unsigned n_node = this->nnode();
630  const unsigned DIM = this->dim();
631 
632  // Find out how many positional dofs there are
633  unsigned n_position_type = this->nnodal_position_type();
634 
635  // Set up memory for the shape/test functions
636  Shape psi(n_node);
637 
638  // Loop over all nodes belonging to the element
639  int local_eqn = 0;
640  // Loop over the nodes of the element
641  for (unsigned l = 0; l < n_node; l++)
642  {
643  // Loop of types of dofs
644  for (unsigned k = 0; k < n_position_type; k++)
645  {
646  // Loop over the force components
647  for (unsigned i = 0; i < DIM; i++)
648  {
649  // Get the local equation
650  local_eqn = this->position_local_eqn(l, k, i);
651  // IF it's not a boundary condition and the interpolated force is nonzero
652  if (local_eqn >= 0)
653  {
654  // add the nodal coupling residual to the global residual vector
655  residuals[local_eqn] += nodal_coupling_residual[l][i];
656  }
657  // Get Jacobian too?
658  if (flag == 1)
659  {
660  // No further action if no coupling forces
661  if (nodal_coupling_jacobian.size() == 0) continue;
662 
663  // Loop over the nodes of the element again
664  for (unsigned ll = 0; ll < n_node; ll++)
665  {
666  // Loop of types of dofs again
667  for (unsigned kk = 0; kk < n_position_type; kk++)
668  {
669  // Loop over the force components again
670  for (unsigned ii = 0; ii < DIM; ii++)
671  {
672  // Get the number of the unknown
673  int local_unknown = this->position_local_eqn(ll, kk, ii);
674  // IF it's not a boundary condition
675  if (( local_unknown >= 0 ) && ( i == ii ) && ( local_unknown >= local_eqn ))
676  {
677  // add the nodal coupling jacobian to the global jacobian matrix
678  jacobian(local_eqn, local_unknown) += nodal_coupling_jacobian[l][ll];
679  if (local_eqn != local_unknown)
680  {
681  jacobian(local_unknown, local_eqn) += nodal_coupling_jacobian[l][ll];
682  }
683  }
684  }
685  }
686  }
687  }
688  }
689  }
690  }
691  }
Vector< Vector< double > > nodal_coupling_residual
Nodal coupling forces.
Definition: VCoupledElement.h:757
Vector< Vector< double > > nodal_coupling_jacobian
Nodal coupling jacobian.
Definition: VCoupledElement.h:760
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

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

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

◆ fill_in_contribution_to_jacobian()

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

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

74  {
75  //Call the generic routine with the flag set to 1
76  fill_in_generic_contribution_to_residuals_pvd(residuals, jacobian, 1);
77 
78  // Add nodal coupling force to the residuals
79  add_internal_coupling_forces_to_residuals(residuals, jacobian, 1);
80  }
void add_internal_coupling_forces_to_residuals(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add the point source contribution to the residual vector.
Definition: VCoupledElement.h:621
void fill_in_generic_contribution_to_residuals_pvd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag) override
Definition: VCoupledElement.h:238

References oomph::VolumeCoupledElement< ELEMENT >::add_internal_coupling_forces_to_residuals(), and oomph::VolumeCoupledElement< ELEMENT >::fill_in_generic_contribution_to_residuals_pvd().

◆ fill_in_contribution_to_residuals()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlineoverride

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

64  {
65  // Call the generic residuals function with flag set to 0 using a dummy matrix argument
66  fill_in_generic_contribution_to_residuals_pvd(residuals, GeneralisedElement::Dummy_matrix, 0);
67 
68  // Add nodal coupling force to the residuals
69  add_internal_coupling_forces_to_residuals(residuals, GeneralisedElement::Dummy_matrix, 0);
70  }

References oomph::VolumeCoupledElement< ELEMENT >::add_internal_coupling_forces_to_residuals(), and oomph::VolumeCoupledElement< ELEMENT >::fill_in_generic_contribution_to_residuals_pvd().

◆ fill_in_generic_contribution_to_residuals_pvd()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::fill_in_generic_contribution_to_residuals_pvd ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
const unsigned &  flag 
)
inlineoverrideprivate
241  {
242 
243  // Get problem dimension
244  const unsigned DIM = this->dim();
245 
246  // Simply set up initial condition?
247  if (this->Solid_ic_pt != 0)
248  {
249  this->fill_in_residuals_for_solid_ic(residuals);
250  return;
251  }
252 
253  //Find out how many nodes there are
254  const unsigned n_node = this->nnode();
255 
256  //Find out how many positional dofs there are
257  const unsigned n_position_type = this->nnodal_position_type();
258 
259  //Set up memory for the shape functions
260  Shape psi(n_node, n_position_type);
261  DShape dpsidxi(n_node, n_position_type, DIM);
262 
263  //Set the value of Nintpt -- the number of integration points
264  const unsigned n_intpt = this->integral_pt()->nweight();
265 
266  //Set the vector to hold the local coordinates in the element
267  Vector<double> s(DIM);
268 
269  // Timescale ratio (non-dim density)
270  double lambda_sq = this->lambda_sq();
271 
272  // Time factor
273  double time_factor = 0.0;
274  if (lambda_sq > 0)
275  {
276  time_factor = this->node_pt(0)->position_time_stepper_pt()->weight(2, 0);
277  }
278 
279  //Integer to store the local equation number
280  int local_eqn = 0;
281 
282  //Loop over the integration points
283  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
284  {
285  //Assign the values of s
286  for (unsigned i = 0; i < DIM; ++i)
287  { s[i] = this->integral_pt()->knot(ipt, i); }
288 
289  //Get the integral weight
290  double w = this->integral_pt()->weight(ipt);
291 
292  //Call the derivatives of the shape functions (and get Jacobian)
293  double J = this->dshape_lagrangian_at_knot(ipt, psi, dpsidxi);
294 
295  //Get the coupling weight at the integration point
296  double coupling_w = 0;
297  for (unsigned l = 0; l < n_node; l++)
298  {
299  double nodal_coupling_w = dynamic_cast<CoupledSolidNode*>(this->node_pt(l))->get_coupling_weight();
300  double psi_ = psi(l);
301  coupling_w += psi_ * nodal_coupling_w;
302  }
303 
304  //Calculate interpolated values of the derivative of global position
305  //wrt lagrangian coordinates
306  DenseMatrix<double> interpolated_G(DIM);
307 
308  // Setup memory for accelerations
309  Vector<double> accel(DIM);
310 
311  //Initialise to zero
312  for (unsigned i = 0; i < DIM; i++)
313  {
314  // Initialise acclerations
315  accel[i] = 0.0;
316  for (unsigned j = 0; j < DIM; j++)
317  {
318  interpolated_G(i, j) = 0.0;
319  }
320  }
321 
322  //Storage for Lagrangian coordinates (initialised to zero)
323  Vector<double> interpolated_xi(DIM, 0.0);
324 
325  //Calculate displacements and derivatives and lagrangian coordinates
326  for (unsigned l = 0; l < n_node; l++)
327  {
328  //Loop over positional dofs
329  for (unsigned k = 0; k < n_position_type; k++)
330  {
331  double psi_ = psi(l, k);
332  //Loop over displacement components (deformed position)
333  for (unsigned i = 0; i < DIM; i++)
334  {
335  //Calculate the Lagrangian coordinates and the accelerations
336  interpolated_xi[i] += this->lagrangian_position_gen(l, k, i) * psi_;
337 
338  // Only compute accelerations if inertia is switched on
339  if (( lambda_sq > 0.0 ) && ( this->Unsteady ))
340  {
341  accel[i] += this->dnodal_position_gen_dt(2, l, k, i) * psi_;
342  }
343 
344  //Loop over derivative directions
345  for (unsigned j = 0; j < DIM; j++)
346  {
347  interpolated_G(j, i) +=
348  this->nodal_position_gen(l, k, i) * dpsidxi(l, k, j);
349  }
350  }
351  }
352  }
353 
354  //Get isotropic growth factor
355  double gamma = 1.0;
356  this->get_isotropic_growth(ipt, s, interpolated_xi, gamma);
357 
358 
359  //Get body force at current time
360  Vector<double> b(DIM);
361  this->body_force(interpolated_xi, b);
362 
363  // We use Cartesian coordinates as the reference coordinate
364  // system. In this case the undeformed metric tensor is always
365  // the identity matrix -- stretched by the isotropic growth
366  double diag_entry = pow(gamma, 2.0 / double(DIM));
367  DenseMatrix<double> g(DIM);
368  for (unsigned i = 0; i < DIM; i++)
369  {
370  for (unsigned j = 0; j < DIM; j++)
371  {
372  if (i == j)
373  { g(i, j) = diag_entry; }
374  else
375  { g(i, j) = 0.0; }
376  }
377  }
378 
379  //Premultiply the undeformed volume ratio (from the isotropic
380  // growth), the integral weights, the coupling weights, and the Jacobian
381  double W = gamma * w * ( 1.0 - coupling_w ) * J;
382 
383  //Declare and calculate the deformed metric tensor
384  DenseMatrix<double> G(DIM);
385 
386  //Assign values of G
387  for (unsigned i = 0; i < DIM; i++)
388  {
389  //Do upper half of matrix
390  for (unsigned j = i; j < DIM; j++)
391  {
392  //Initialise G(i,j) to zero
393  G(i, j) = 0.0;
394  //Now calculate the dot product
395  for (unsigned k = 0; k < DIM; k++)
396  {
397  G(i, j) += interpolated_G(i, k) * interpolated_G(j, k);
398  }
399  }
400  //Matrix is symmetric so just copy lower half
401  for (unsigned j = 0; j < i; j++)
402  {
403  G(i, j) = G(j, i);
404  }
405  }
406 
407  //Now calculate the stress tensor from the constitutive law
408  DenseMatrix<double> sigma(DIM);
409  ELEMENT::get_stress(g, G, sigma);
410 
411  // Add pre-stress
412  for (unsigned i = 0; i < DIM; i++)
413  {
414  for (unsigned j = 0; j < DIM; j++)
415  {
416  sigma(i, j) += this->prestress(i, j, interpolated_xi);
417  }
418  }
419 
420  // Get stress derivative by FD only needed for Jacobian
421  //-----------------------------------------------------
422 
423  // Stress derivative
424  RankFourTensor<double> d_stress_dG(DIM, DIM, DIM, DIM, 0.0);
425  // Derivative of metric tensor w.r.t. to nodal coords
426  RankFiveTensor<double> d_G_dX(n_node, n_position_type, DIM, DIM, DIM, 0.0);
427 
428  // Get Jacobian too?
429  if (flag == 1)
430  {
431  // Derivative of metric tensor w.r.t. to discrete positional dofs
432  // NOTE: Since G is symmetric we only compute the upper triangle
433  // and DO NOT copy the entries across. Subsequent computations
434  // must (and, in fact, do) therefore only operate with upper
435  // triangular entries
436  for (unsigned ll = 0; ll < n_node; ll++)
437  {
438  for (unsigned kk = 0; kk < n_position_type; kk++)
439  {
440  for (unsigned ii = 0; ii < DIM; ii++)
441  {
442  for (unsigned aa = 0; aa < DIM; aa++)
443  {
444  for (unsigned bb = aa; bb < DIM; bb++)
445  {
446  d_G_dX(ll, kk, ii, aa, bb) =
447  interpolated_G(aa, ii) * dpsidxi(ll, kk, bb) +
448  interpolated_G(bb, ii) * dpsidxi(ll, kk, aa);
449  }
450  }
451  }
452  }
453  }
454 
455  //Get the "upper triangular" entries of the derivatives of the stress
456  //tensor with respect to G
457  this->get_d_stress_dG_upper(g, G, sigma, d_stress_dG);
458  }
459 
460  //=====EQUATIONS OF ELASTICITY FROM PRINCIPLE OF VIRTUAL DISPLACEMENTS========
461 
462  //Loop over the test functions, nodes of the element
463  for (unsigned l = 0; l < n_node; l++)
464  {
465  //Loop of types of dofs
466  for (unsigned k = 0; k < n_position_type; k++)
467  {
468  // Offset for faster access
469  const unsigned offset5 = dpsidxi.offset(l, k);
470 
471  //Loop over the displacement components
472  for (unsigned i = 0; i < DIM; i++)
473  {
474  //Get the equation number
475  local_eqn = this->position_local_eqn(l, k, i);
476 
477  /*IF it's not a boundary condition*/
478  if (local_eqn >= 0)
479  {
480  //Initialise contribution to sum
481  double sum = 0.0;
482 
483  // Acceleration and body force
484  sum += ( lambda_sq * accel[i] - b[i] ) * psi(l, k);
485 
486  // Stress term
487  for (unsigned a = 0; a < DIM; a++)
488  {
489  unsigned count = offset5;
490  for (unsigned b = 0; b < DIM; b++)
491  {
492  //Add the stress terms to the residuals
493  sum += sigma(a, b) * interpolated_G(a, i) *
494  dpsidxi.raw_direct_access(count);
495  ++count;
496  }
497  }
498  residuals[local_eqn] += W * sum;
499 
500  // Get Jacobian too?
501  if (flag == 1)
502  {
503  // Offset for faster access in general stress loop
504  const unsigned offset1 = d_G_dX.offset(l, k, i);
505 
506  //Loop over the nodes of the element again
507  for (unsigned ll = 0; ll < n_node; ll++)
508  {
509  //Loop of types of dofs again
510  for (unsigned kk = 0; kk < n_position_type; kk++)
511  {
512  //Loop over the displacement components again
513  for (unsigned ii = 0; ii < DIM; ii++)
514  {
515  //Get the number of the unknown
516  int local_unknown = this->position_local_eqn(ll, kk, ii);
517 
518  /*IF it's not a boundary condition*/
519  if (local_unknown >= 0)
520  {
521  // Offset for faster access in general stress loop
522  const unsigned offset2 = d_G_dX.offset(ll, kk, ii);
523  const unsigned offset4 = dpsidxi.offset(ll, kk);
524 
525  // General stress term
526  //--------------------
527  double sum = 0.0;
528  unsigned count1 = offset1;
529  for (unsigned a = 0; a < DIM; a++)
530  {
531  // Bump up direct access because we're only
532  // accessing upper triangle
533  count1 += a;
534  for (unsigned b = a; b < DIM; b++)
535  {
536  double factor = d_G_dX.raw_direct_access(count1);
537  if (a == b) factor *= 0.5;
538 
539  // Offset for faster access
540  unsigned offset3 = d_stress_dG.offset(a, b);
541  unsigned count2 = offset2;
542  unsigned count3 = offset3;
543 
544  for (unsigned aa = 0; aa < DIM; aa++)
545  {
546  // Bump up direct access because we're only
547  // accessing upper triangle
548  count2 += aa;
549  count3 += aa;
550 
551  // Only upper half of derivatives w.r.t. symm tensor
552  for (unsigned bb = aa; bb < DIM; bb++)
553  {
554  sum += factor *
555  d_stress_dG.raw_direct_access(count3) *
556  d_G_dX.raw_direct_access(count2);
557  ++count2;
558  ++count3;
559  }
560  }
561  ++count1;
562  }
563 
564  }
565 
566  // Multiply by weight and add contribution
567  // (Add directly because this bit is nonsymmetric)
568  jacobian(local_eqn, local_unknown) += sum * W;
569 
570  // Only upper triangle (no separate test for bc as
571  // local_eqn is already nonnegative)
572  if (( i == ii ) && ( local_unknown >= local_eqn ))
573  {
574  //Initialise contribution
575  double sum = 0.0;
576 
577  // Inertia term
578  sum += lambda_sq * time_factor * psi(ll, kk) * psi(l, k);
579 
580  // Stress term
581  unsigned count4 = offset4;
582  for (unsigned a = 0; a < DIM; a++)
583  {
584  //Cache term
585  const double factor =
586  dpsidxi.raw_direct_access(count4);// ll ,kk
587  ++count4;
588 
589  unsigned count5 = offset5;
590  for (unsigned b = 0; b < DIM; b++)
591  {
592  sum += sigma(a, b) * factor *
593  dpsidxi.raw_direct_access(count5); // l ,k
594  ++count5;
595  }
596  }
597  //Add contribution to jacobian
598  jacobian(local_eqn, local_unknown) += sum * W;
599  //Add to lower triangular section
600  if (local_eqn != local_unknown)
601  {
602  jacobian(local_unknown, local_eqn) += sum * W;
603  }
604  }
605 
606  } //End of if not boundary condition
607  }
608  }
609  }
610  }
611 
612  } //End of if not boundary condition
613 
614  } //End of loop over coordinate directions
615  } //End of loop over type of dof
616  } //End of loop over shape functions
617  } //End of loop over integration points
618  }
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.

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

◆ get_momentum_and_energy()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::get_momentum_and_energy ( double mass,
Vector< double > &  lin_mo,
Vector< double > &  ang_mo,
double pot_en,
double kin_en 
)
inline
102  {
103  const unsigned DIM = this->dim();
104  // Initialise mass
105  mass = 0;
106  // Initialise momentum
107  lin_mo.initialise(0);
108  ang_mo.initialise(0);
109  // Initialise energy
110  pot_en = 0;
111  kin_en = 0;
112 
113  //Set the value of n_intpt
114  unsigned n_intpt = this->integral_pt()->nweight();
115 
116  //Set the Vector to hold local coordinates
117  Vector<double> s(DIM);
118 
119  //Find out how many nodes there are
120  const unsigned n_node = this->nnode();
121 
122  //Find out how many positional dofs there are
123  const unsigned n_position_type = this->nnodal_position_type();
124 
125  //Set up memory for the shape functions
126  Shape psi(n_node, n_position_type);
127  DShape dpsidxi(n_node, n_position_type, DIM);
128 
129  // Timescale ratio (non-dim density)
130  double lambda_sq = this->lambda_sq();
131 
132  //Loop over the integration points
133  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
134  {
135  //Assign values of s
136  for (unsigned i = 0; i < DIM; i++)
137  { s[i] = this->integral_pt()->knot(ipt, i); }
138 
139  //Get the integral weight
140  double w = this->integral_pt()->weight(ipt);
141 
142  //Call the derivatives of the shape functions and get Jacobian
143  double J = this->dshape_lagrangian_at_knot(ipt, psi, dpsidxi);
144 
145  //Get mass and the coupling weight at the integration point
146  double coupling_w = 0;
147  for (unsigned l = 0; l < n_node; l++)
148  {
149  double nodal_coupling_w = dynamic_cast<CoupledSolidNode*>(this->node_pt(l))->get_coupling_weight();
150  double psi_ = psi(l);
151  coupling_w += psi_ * nodal_coupling_w;
152  }
153 
154  //Storage for Lagrangian coordinates and velocity (initialised to zero)
155  Vector<double> interpolated_xi(DIM, 0.0);
156  Vector<double> veloc(DIM, 0.0);
157 
158  //Calculate lagrangian coordinates
159  for (unsigned l = 0; l < n_node; l++)
160  {
161  //Loop over positional dofs
162  for (unsigned k = 0; k < n_position_type; k++)
163  {
164  //Loop over displacement components (deformed position)
165  for (unsigned i = 0; i < DIM; i++)
166  {
167  //Calculate the Lagrangian coordinates
168  interpolated_xi[i] += this->lagrangian_position_gen(l, k, i) * psi(l, k);
169 
170  //Calculate the velocity components (if unsteady solve)
171  if (this->Unsteady)
172  {
173  veloc[i] += this->dnodal_position_gen_dt(l, k, i) * psi(l, k);
174  }
175  }
176  }
177  }
178 
179  //Get isotropic growth factor
180  double gamma = 1.0;
181  this->get_isotropic_growth(ipt, s, interpolated_xi, gamma);
182 
183  //Premultiply the undeformed volume ratio (from the isotropic
184  // growth), the integral weights, the coupling weights, and the Jacobian
185  double W = gamma * w * ( 1.0 - coupling_w ) * J;
186 
187  DenseMatrix<double> sigma(DIM, DIM);
188  DenseMatrix<double> strain(DIM, DIM);
189 
190  //Now calculate the stress tensor from the constitutive law
191  this->get_stress(s, sigma);
192 
193  // Add pre-stress
194  for (unsigned i = 0; i < DIM; i++)
195  {
196  for (unsigned j = 0; j < DIM; j++)
197  {
198  sigma(i, j) += this->prestress(i, j, interpolated_xi);
199  }
200  }
201 
202  //get the strain
203  this->get_strain(s, strain);
204 
205  // Initialise
206  double local_pot_en = 0;
207  double veloc_sq = 0;
208 
209  // Compute integrals
210  for (unsigned i = 0; i < DIM; i++)
211  {
212  for (unsigned j = 0; j < DIM; j++)
213  {
214  local_pot_en += sigma(i, j) * strain(i, j);
215  }
216  veloc_sq += veloc[i] * veloc[i];
217  }
218 
219  // Mass
220  mass += lambda_sq * W;
221  // Linear momentum and angular momentum
222  Vector<double> cross_product(DIM, 0);
223  VectorHelpers::cross(interpolated_xi, veloc, cross_product);
224  for (unsigned i = 0; i < DIM; i++)
225  {
226  lin_mo[i] += lambda_sq * veloc[i] * W;
227  ang_mo[i] += lambda_sq * cross_product[i] * W;
228  }
229  // Potential energy
230  pot_en += 0.5 * local_pot_en * W;
231  // Kinetic energy
232  kin_en += lambda_sq * 0.5 * veloc_sq * W;
233  }
234  }

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

◆ get_nodal_coupling_jacobian()

template<class ELEMENT >
double& oomph::VolumeCoupledElement< ELEMENT >::get_nodal_coupling_jacobian ( const unsigned &  l,
const unsigned &  ll 
)
inline

◆ get_nodal_coupling_residual()

template<class ELEMENT >
double& oomph::VolumeCoupledElement< ELEMENT >::get_nodal_coupling_residual ( const unsigned &  l,
const unsigned &  i 
)
inline

◆ getCouplingStiffness()

template<class ELEMENT >
double& oomph::VolumeCoupledElement< ELEMENT >::getCouplingStiffness ( const unsigned &  m,
const unsigned &  l 
)
inline
98  { return couplingStiffness[m][l]; }
Vector< Vector< double > > couplingStiffness
Coupling stiffness to discrete particles (FIXME: should be moved into DPMVCoupledElement)
Definition: VCoupledElement.h:763

References oomph::VolumeCoupledElement< ELEMENT >::couplingStiffness.

◆ output()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::output ( std::ostream &  outfile,
const unsigned &  n_plot 
)
inlineoverrideprivate
694  {
695  const unsigned DIM = this->dim();
696  Vector<double> x(DIM);
697  Vector<double> dxdt(DIM);
698  Vector<double> s(DIM);
699 
700  // Tecplot header info
701  outfile << this->tecplot_zone_string(n_plot);
702 
703  // Loop over plot points
704  unsigned num_plot_points = this->nplot_points(n_plot);
705  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
706  {
707  // Get local coordinates of plot point
708  this->get_s_plot(iplot, n_plot, s);
709 
710  // Get Eulerian coordinates
711  this->interpolated_x(s, x);
712  SolidFiniteElement* el_pt = dynamic_cast<SolidFiniteElement*>(this);
713  el_pt->interpolated_dxdt(s, 1, dxdt);
714 
715  // Dummy integration point
716  unsigned ipt = 0;
717 
718  // Output the x,y,..
719  for (unsigned i = 0; i < DIM; i++)
720  { outfile << x[i] * Global_Physical_Variables::lenScale << " "; }
721 
722  // Output velocity dnodal_position_gen_dt
723  for (unsigned i = 0; i < DIM; i++)
724  {
726  << " ";
727  }
728 
729  //Find the number of nodes
730  const unsigned n_node = this->nnode();
731  //Find the number of positional types
732  const unsigned n_position_type = this->nnodal_position_type();
733  //Assign storage for the local shape function
734  Shape psi(n_node, n_position_type);
735  //Find the values of shape function
736  this->shape(s, psi);
737  // Initialize coupling weight
738  double w = 0;
739  // Loop over the local nodes
740  for (unsigned l = 0; l < n_node; l++)
741  {
742  double nodal_coupling_w = dynamic_cast<CoupledSolidNode*>(this->node_pt(l))->get_coupling_weight();
743  w += nodal_coupling_w * psi(l);
744  }
745  // Output coupling weight
746  outfile << w << " ";
747 
748  outfile << std::endl;
749  }
750 
751  // Write tecplot footer (e.g. FE connectivity lists)
752  this->write_tecplot_zone_footer(outfile, n_plot);
753  outfile << std::endl;
754  }
double lenScale
Length scale.
Definition: VCoupledElement.h:31
double timeScale
Time scale.
Definition: VCoupledElement.h:34

References constants::i, Global_Physical_Variables::lenScale, and Global_Physical_Variables::timeScale.

◆ set_nodal_coupling_jacobian()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::set_nodal_coupling_jacobian ( Vector< Vector< double > > &  cJacobian)
inline

◆ set_nodal_coupling_residual()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::set_nodal_coupling_residual ( Vector< Vector< double > > &  cResidual)
inline

◆ setCouplingStiffness()

template<class ELEMENT >
void oomph::VolumeCoupledElement< ELEMENT >::setCouplingStiffness ( Vector< Vector< double > > &  cStiffness)
inline

Member Data Documentation

◆ couplingStiffness

template<class ELEMENT >
Vector<Vector<double> > oomph::VolumeCoupledElement< ELEMENT >::couplingStiffness
private

Coupling stiffness to discrete particles (FIXME: should be moved into DPMVCoupledElement)

Referenced by oomph::VolumeCoupledElement< ELEMENT >::getCouplingStiffness(), and oomph::VolumeCoupledElement< ELEMENT >::setCouplingStiffness().

◆ nodal_coupling_jacobian

◆ nodal_coupling_residual


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