refineable_solid_elements.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2022 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Header file for refineable solid mechanics elements
27 
28 // Include guard to prevent multiple inclusions of this header
29 #ifndef OOMPH_REFINEABLE_ELASTICITY_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_ELASTICITY_ELEMENTS_HEADER
31 
32 // oomph-lib headers
33 #include "solid_elements.h"
34 #include "../generic/refineable_quad_element.h"
35 #include "../generic/refineable_brick_element.h"
36 #include "../generic/error_estimator.h"
37 
38 namespace oomph
39 {
40  //========================================================================
42  //========================================================================
43  template<unsigned DIM>
44  class RefineablePVDEquations : public virtual PVDEquations<DIM>,
45  public virtual RefineableSolidElement,
46  public virtual ElementWithZ2ErrorEstimator
47  {
48  public:
51  : PVDEquations<DIM>(),
55  {
56  }
57 
60  Vector<double>& residuals,
61  DenseMatrix<double>& jacobian,
62  const unsigned& flag);
63 
65  void get_interpolated_values(const unsigned& t,
66  const Vector<double>& s,
67  Vector<double>& values)
68  {
69  values.clear();
70  }
71 
74  Vector<double>& values)
75  {
76  values.clear();
77  }
78 
80  unsigned num_Z2_flux_terms()
81  {
82  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
83  return DIM + DIM * (DIM - 1) / 2;
84  }
85 
89  {
90 #ifdef PARANOID
91  unsigned num_entries = DIM + ((DIM * DIM) - DIM) / 2;
92  if (flux.size() != num_entries)
93  {
94  std::ostringstream error_message;
95  error_message << "The flux vector has the wrong number of entries, "
96  << flux.size() << ", whereas it should be " << num_entries
97  << std::endl;
98  throw OomphLibError(error_message.str(),
101  }
102 #endif
103 
104  // Get strain matrix
105  DenseMatrix<double> strain(DIM);
106  this->get_strain(s, strain);
107 
108  // Pack into flux Vector
109  unsigned icount = 0;
110 
111  // Start with diagonal terms
112  for (unsigned i = 0; i < DIM; i++)
113  {
114  flux[icount] = strain(i, i);
115  icount++;
116  }
117 
118  // Off diagonals row by row
119  for (unsigned i = 0; i < DIM; i++)
120  {
121  for (unsigned j = i + 1; j < DIM; j++)
122  {
123  flux[icount] = strain(i, j);
124  icount++;
125  }
126  }
127  }
128 
130  unsigned ncont_interpolated_values() const
131  {
132  return 0;
133  }
134 
135  // Return a pointer to the solid node at which pressure dof l2 is stored
136  // This is only required so that the generic templating in
137  // PseudoSolidNodeUpdateElements works OK
138  virtual Node* solid_pressure_node_pt(const unsigned& l)
139  {
140  return 0;
141  }
142 
145  {
146  RefineablePVDEquations<DIM>* cast_father_element_pt =
147  dynamic_cast<RefineablePVDEquations<DIM>*>(this->father_element_pt());
148 
149  // Do whatever needs to be done in the base class
151 
152  // Set pointer to isotropic growth function
154  cast_father_element_pt->isotropic_growth_fct_pt();
155 
156  // Set pointer to body force function
157  this->Body_force_fct_pt = cast_father_element_pt->body_force_fct_pt();
158 
159  // Set pointer to the contitutive law
160  this->Constitutive_law_pt = cast_father_element_pt->constitutive_law_pt();
161 
162  // Set the timescale ratio (non-dim. density)
163  this->Lambda_sq_pt = cast_father_element_pt->lambda_sq_pt();
164 
165  // Set the flag that switches inertia on/off
166  this->Unsteady = cast_father_element_pt->is_inertia_enabled();
167 
168  // Evaluation of Jacobian by same method as father
170  cast_father_element_pt->is_jacobian_evaluated_by_fd();
171  }
172  };
173 
174  //========================================================================
176  //========================================================================
177  template<unsigned DIM, unsigned NNODE_1D>
178  class RefineableQPVDElement : public virtual QPVDElement<DIM, NNODE_1D>,
179  public virtual RefineablePVDEquations<DIM>,
180  public virtual RefineableSolidQElement<DIM>
181  {
182  public:
185  : QPVDElement<DIM, NNODE_1D>(),
190  {
191  }
192 
194  void rebuild_from_sons(Mesh*& mesh_pt) {}
195 
197  unsigned nvertex_node() const
198  {
200  }
201 
203  Node* vertex_node_pt(const unsigned& j) const
204  {
206  }
207 
210  unsigned nrecovery_order()
211  {
212  return NNODE_1D - 1;
213  }
214 
218  };
219 
220  //==============================================================
222  //==============================================================
223  template<unsigned NNODE_1D>
225  : public virtual SolidQElement<1, NNODE_1D>
226  {
227  public:
228  // Make sure that we call the constructor of the SolidQElement
229  // Only the Intel compiler seems to need this!
230  FaceGeometry() : SolidQElement<1, NNODE_1D>() {}
231  };
232 
233  //==============================================================
235  //==============================================================
236  template<unsigned NNODE_1D>
238  : public virtual PointElement
239  {
240  public:
241  // Make sure that we call the constructor of the SolidQElement
242  // Only the Intel compiler seems to need this!
244  };
245 
246 
247  //==============================================================
249  //==============================================================
250  template<unsigned NNODE_1D>
252  : public virtual SolidQElement<2, NNODE_1D>
253  {
254  public:
255  // Make sure that we call the constructor of the SolidQElement
256  // Only the Intel compiler seems to need this!
257  FaceGeometry() : SolidQElement<2, NNODE_1D>() {}
258  };
259 
260  //==============================================================
262  //==============================================================
263  template<unsigned NNODE_1D>
265  : public virtual SolidQElement<1, NNODE_1D>
266  {
267  public:
268  // Make sure that we call the constructor of the SolidQElement
269  // Only the Intel compiler seems to need this!
270  FaceGeometry() : SolidQElement<1, NNODE_1D>() {}
271  };
272 
273 
274  //===========================================================================
278  //===========================================================================
279  template<unsigned DIM>
281  : public virtual PVDEquationsWithPressure<DIM>,
282  public virtual RefineableSolidElement,
283  public virtual ElementWithZ2ErrorEstimator
284  {
285  public:
292  {
293  }
294 
300  Vector<double>& residuals,
301  DenseMatrix<double>& jacobian,
302  DenseMatrix<double>& mass_matrix,
303  const unsigned& flag);
304 
306  void get_interpolated_values(const unsigned& t,
307  const Vector<double>& s,
308  Vector<double>& values)
309  {
310  values.clear();
311  }
312 
315  Vector<double>& values)
316  {
317  values.clear();
318  }
319 
321  unsigned num_Z2_flux_terms()
322  {
323  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
324  return DIM + DIM * (DIM - 1) / 2;
325  }
326 
327  // Get 'flux' for Z2 error recovery: Upper triangular entries
329  //----------------------------------------------------------------
331  {
332  // Find the dimension of the problem
333 #ifdef PARANOID
334  unsigned num_entries = DIM + ((DIM * DIM) - DIM) / 2;
335  if (flux.size() != num_entries)
336  {
337  std::ostringstream error_message;
338  error_message << "The flux vector has the wrong number of entries, "
339  << flux.size() << ", whereas it should be " << num_entries
340  << std::endl;
341  throw OomphLibError(error_message.str(),
344  }
345 #endif
346 
347  // Get strain matrix
348  DenseMatrix<double> strain(DIM);
349  this->get_strain(s, strain);
350 
351  // Pack into flux Vector
352  unsigned icount = 0;
353 
354  // Start with diagonal terms
355  for (unsigned i = 0; i < DIM; i++)
356  {
357  flux[icount] = strain(i, i);
358  icount++;
359  }
360 
361  // Off diagonals row by row
362  for (unsigned i = 0; i < DIM; i++)
363  {
364  for (unsigned j = i + 1; j < DIM; j++)
365  {
366  flux[icount] = strain(i, j);
367  icount++;
368  }
369  }
370  }
371 
373  unsigned ncont_interpolated_values() const
374  {
375  return 0;
376  }
377 
378  // Return a pointer to the solid node at which pressure dof l2 is stored
379  virtual Node* solid_pressure_node_pt(const unsigned& l)
380  {
381  return 0;
382  }
383 
384 
387  {
388  RefineablePVDEquationsWithPressure<DIM>* cast_father_element_pt =
390  this->father_element_pt());
391 
392  // Do whatever needs to be done in the base class
394 
395  // Set pointer to isotropic growth function
397  cast_father_element_pt->isotropic_growth_fct_pt();
398 
399  // Set pointer to body force function
400  this->Body_force_fct_pt = cast_father_element_pt->body_force_fct_pt();
401 
402  // Set pointer to the contitutive law
403  this->Constitutive_law_pt = cast_father_element_pt->constitutive_law_pt();
404 
405  // Set the timescale ratio (non-dim. density)
406  this->Lambda_sq_pt = cast_father_element_pt->lambda_sq_pt();
407 
408  // Set the flag that switches inertia on/off
409  this->Unsteady = cast_father_element_pt->is_inertia_enabled();
410 
411  // Set the incompressibility flag
412  this->Incompressible = cast_father_element_pt->is_incompressible();
413 
414  // Evaluation of Jacobian by same method as father
416  cast_father_element_pt->is_jacobian_evaluated_by_fd();
417  }
418 
419 
422  void get_mass_matrix_diagonal(Vector<double>& mass_diag);
423  };
424 
425  //===========================================================================
430  //===========================================================================
431  template<unsigned DIM>
433  : public virtual QPVDElementWithPressure<DIM>,
434  public virtual RefineablePVDEquationsWithPressure<DIM>,
435  public virtual RefineableSolidQElement<DIM>
436  {
437  private:
440  {
441  unsigned n_pres = this->npres_solid();
442  // loop over pressure dofs and unpin them
443  for (unsigned l = 0; l < n_pres; l++)
444  {
446  }
447  }
448 
449  public:
457  {
458  }
459 
460 
463  inline void rebuild_from_sons(Mesh*& mesh_pt);
464 
466  unsigned nvertex_node() const
467  {
469  }
470 
472  Node* vertex_node_pt(const unsigned& j) const
473  {
475  }
476 
479  unsigned nrecovery_order()
480  {
481  return 2;
482  } // NNODE_1D-1;}
483 
487 
490  inline void further_build();
491 
492 
494  unsigned ncont_interpolated_values() const
495  {
496  return 0;
497  }
498  };
499 
500  //======================================================================
502  //=======================================================================
503  template<>
505  : public virtual SolidQElement<1, 3>
506  {
507  public:
508  // Make sure that we call the constructor of the SolidQElement
509  // Only the Intel compiler seems to need this!
511  };
512 
513 
514  //===========================================================================
517  //============================================================================
518  template<>
520  : public virtual PointElement
521  {
522  public:
523  // Make sure that we call the constructor of the SolidQElement
524  // Only the Intel compiler seems to need this!
526  };
527 
528  //====================================================================
530  //===================================================================
531  template<>
533  Mesh*& mesh_pt)
534  {
535  using namespace QuadTreeNames;
536 
537  // Storage for the solid pressure of each son
538  double centre_solid_p[4] = {0.0, 0.0, 0.0, 0.0};
539 
540  // Loop over the sons and assign the central solid pressure
541  for (unsigned ison = 0; ison < 4; ison++)
542  {
543  // Add the sons midnode pressure
544  centre_solid_p[ison] =
546  quadtree_pt()->son_pt(ison)->object_pt())
547  ->solid_p(0);
548  }
549 
550  // Use the average for the central solid pressure
551  double p_value = 0.25 * (centre_solid_p[0] + centre_solid_p[1] +
552  centre_solid_p[2] + centre_solid_p[3]);
553  // Actually set the pressure
554  set_solid_p(0, p_value);
555 
556 
557  // Slope in s_0 direction
558  //----------------------
559 
560  // Use average of the 2 FD approximations based on the
561  // elements central pressure values
562  // [Other options: Take average of the four
563  // pressure derivatives]
564  double slope1 = centre_solid_p[SE] - centre_solid_p[SW];
565  double slope2 = centre_solid_p[NE] - centre_solid_p[NW];
566 
567 
568  // Use the average value
569  p_value = 0.5 * (slope1 + slope2);
570  set_solid_p(1, p_value);
571 
572  // Slope in s_1 direction
573  //----------------------
574 
575  // Use average of the 2 FD approximations based on the
576  // elements central pressure values
577  // [Other options: Take average of the four
578  // pressure derivatives]
579 
580  slope1 = centre_solid_p[NE] - centre_solid_p[SE];
581  slope2 = centre_solid_p[NW] - centre_solid_p[SW];
582 
583  // Use the average
584  p_value = 0.5 * (slope1 + slope2);
585  set_solid_p(2, p_value);
586  }
587 
588 
589  //==============================================================
592  //=============================================================
593  template<>
595  {
597 
598  using namespace QuadTreeNames;
599 
600  // What type of son am I? Ask my quadtree representation...
601  int son_type = quadtree_pt()->son_type();
602 
603  // Pointer to my father (in element impersonation)
604  RefineableElement* father_el_pt = quadtree_pt()->father_pt()->object_pt();
605 
606  Vector<double> s_father(2);
607 
608  // Son midpoint is located at the following coordinates in father element:
609 
610  // South west son
611  if (son_type == SW)
612  {
613  s_father[0] = -0.5;
614  s_father[1] = -0.5;
615  }
616  // South east son
617  else if (son_type == SE)
618  {
619  s_father[0] = 0.5;
620  s_father[1] = -0.5;
621  }
622  // North east son
623  else if (son_type == NE)
624  {
625  s_father[0] = 0.5;
626  s_father[1] = 0.5;
627  }
628 
629  // North west son
630  else if (son_type == NW)
631  {
632  s_father[0] = -0.5;
633  s_father[1] = 0.5;
634  }
635 
636  // Pressure value in father element
637  RefineableQPVDElementWithPressure<2>* cast_father_element_pt =
638  dynamic_cast<RefineableQPVDElementWithPressure<2>*>(father_el_pt);
639  double press = cast_father_element_pt->interpolated_solid_p(s_father);
640 
641  // Pressure value gets copied straight into internal dof:
642  set_solid_p(0, press);
643 
644  // The slopes get copied from father and halved
645  set_solid_p(1, 0.5 * cast_father_element_pt->solid_p(1));
646  set_solid_p(2, 0.5 * cast_father_element_pt->solid_p(2));
647  }
648 
649 
650  //===============================================================
652  //===============================================================
653  template<>
655  Mesh*& mesh_pt)
656  {
657  using namespace OcTreeNames;
658 
659  // Storage for the central solid pressure of each son
660  double centre_solid_p[8] = {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0};
661 
662  // Loop over the sons and assign the central solid pressure
663  for (unsigned ison = 0; ison < 8; ison++)
664  {
665  centre_solid_p[ison] = octree_pt()
666  ->son_pt(ison)
667  ->object_pt()
668  ->internal_data_pt(this->P_solid_internal_index)
669  ->value(0);
670  }
671 
672  // Central pressure value:
673  //-----------------------
674 
675  // Use average of the sons central pressure values
676  // Other options: Take average of the four (discontinuous)
677  // pressure values at the father's midpoint]
678  double av_press = 0.0;
679 
680  // Loop over the sons and sum the centre pressures
681  for (unsigned ison = 0; ison < 8; ison++)
682  {
683  av_press += centre_solid_p[ison];
684  }
685 
686  // Use the average
687  internal_data_pt(this->P_solid_internal_index)
688  ->set_value(0, 0.125 * av_press);
689 
690 
691  // Slope in s_0 direction
692  //----------------------
693 
694  // Use average of the 4 FD approximations based on the
695  // elements central pressure values
696  // [Other options: Take average of the four
697  // pressure derivatives]
698 
699  double slope1 = centre_solid_p[RDF] - centre_solid_p[LDF];
700  double slope2 = centre_solid_p[RUF] - centre_solid_p[LUF];
701  double slope3 = centre_solid_p[RDB] - centre_solid_p[LDB];
702  double slope4 = centre_solid_p[RUB] - centre_solid_p[LUB];
703 
704  // Use the average
705  internal_data_pt(this->P_solid_internal_index)
706  ->set_value(1, 0.25 * (slope1 + slope2 + slope3 + slope4));
707 
708 
709  // Slope in s_1 direction
710  //----------------------
711 
712  // Use average of the 4 FD approximations based on the
713  // elements central pressure values
714  // [Other options: Take average of the four
715  // pressure derivatives]
716  slope1 = centre_solid_p[LUB] - centre_solid_p[LDB];
717  slope2 = centre_solid_p[RUB] - centre_solid_p[RDB];
718  slope3 = centre_solid_p[LUF] - centre_solid_p[LDF];
719  slope4 = centre_solid_p[RUF] - centre_solid_p[RDF];
720 
721  // Use the average
722  internal_data_pt(this->P_solid_internal_index)
723  ->set_value(2, 0.25 * (slope1 + slope2 + slope3 + slope4));
724 
725 
726  // Slope in s_2 direction
727  //----------------------
728 
729  // Use average of the 4 FD approximations based on the
730  // elements central pressure values
731  // [Other options: Take average of the four
732  // pressure derivatives]
733  slope1 = centre_solid_p[LUF] - centre_solid_p[LUB];
734  slope2 = centre_solid_p[RUF] - centre_solid_p[RUB];
735  slope3 = centre_solid_p[LDF] - centre_solid_p[LDB];
736  slope4 = centre_solid_p[RDF] - centre_solid_p[RDB];
737 
738  // Use the average
739  internal_data_pt(this->P_solid_internal_index)
740  ->set_value(3, 0.25 * (slope1 + slope2 + slope3 + slope4));
741  }
742 
743 
744  //================================================================
746  //=================================================================
747  template<>
749 
750  {
752 
753  using namespace OcTreeNames;
754 
755  // What type of son am I? Ask my quadtree representation...
756  int son_type = octree_pt()->son_type();
757 
758  // Pointer to my father (in element impersonation)
759  RefineableQElement<3>* father_el_pt = dynamic_cast<RefineableQElement<3>*>(
760  octree_pt()->father_pt()->object_pt());
761 
762  Vector<double> s_father(3);
763 
764  // Son midpoint is located at the following coordinates in father element:
765  for (unsigned i = 0; i < 3; i++)
766  {
767  s_father[i] = 0.5 * OcTree::Direction_to_vector[son_type][i];
768  }
769 
770  // Pressure value in father element
771  RefineableQPVDElementWithPressure<3>* cast_father_element_pt =
772  dynamic_cast<RefineableQPVDElementWithPressure<3>*>(father_el_pt);
773 
774  double press = cast_father_element_pt->interpolated_solid_p(s_father);
775 
776 
777  // Pressure value gets copied straight into internal dof:
778  set_solid_p(0, press);
779 
780  // The slopes get copied from father and halved
781  for (unsigned i = 1; i < 4; i++)
782  {
783  set_solid_p(i, 0.5 * cast_father_element_pt->solid_p(i));
784  }
785  }
786 
787  //=========================================================================
789  //========================================================================
790  template<>
792  : public virtual SolidQElement<2, 3>
793  {
794  public:
795  // Make sure that we call the constructor of the SolidQElement
796  // Only the Intel compiler seems to need this!
798  };
799 
800 
801  //========================================================================
804  //==========================================================================
805  template<>
807  : public virtual SolidQElement<1, 3>
808  {
809  public:
810  // Make sure that we call the constructor of the SolidQElement
811  // Only the Intel compiler seems to need this!
813  };
814 
815 
816  //===========================================================================
820  //===========================================================================
821  template<unsigned DIM>
823  : public virtual QPVDElementWithContinuousPressure<DIM>,
824  public virtual RefineablePVDEquationsWithPressure<DIM>,
825  public virtual RefineableSolidQElement<DIM>
826  {
827  public:
835  {
836  }
837 
838 
841  unsigned required_nvalue(const unsigned& n) const
842  {
843  return 1;
844  }
845 
847  unsigned ncont_interpolated_values() const
848  {
849  return 1;
850  }
851 
853  void rebuild_from_sons(Mesh*& mesh_pt) {}
854 
857  Vector<double>& values)
858  {
859  // There is only one solid pressure, initialise to zero
860  values.resize(1);
861 
862  // Get the interpolated value
863  values[0] = this->interpolated_solid_p(s);
864  }
865 
867  void get_interpolated_values(const unsigned& t,
868  const Vector<double>& s,
869  Vector<double>& values)
870  {
871  // There is only one solid pressure, initialise to zero
872  values.resize(1);
873  // The solid pressure does not depend on time!
874  values[0] = this->interpolated_solid_p(s);
875  }
876 
877 
880  {
881  // find the index at which the pressure is stored
882  int solid_p_index = this->solid_p_nodal_index();
883  unsigned n_node = this->nnode();
884  // loop over nodes
885  for (unsigned n = 0; n < n_node; n++)
886  {
887  this->node_pt(n)->unpin(solid_p_index);
888  }
889  }
890 
891 
894  {
895  // Find the index of the solid pressure
896  int solid_p_index = this->solid_p_nodal_index();
897  // Let's pin all pressure nodes
898  unsigned n_node = this->nnode();
899  for (unsigned l = 0; l < n_node; l++)
900  {
901  // Pin the solid pressure
902  this->node_pt(l)->pin(solid_p_index);
903  }
904 
905  // Now loop over the pressure nodes and unpin the solid pressures
906  unsigned n_solid_pres = this->npres_solid();
907  // Loop over these nodes and unpin the solid pressures
908  for (unsigned l = 0; l < n_solid_pres; l++)
909  {
910  Node* nod_pt = this->solid_pressure_node_pt(l);
911  if (!nod_pt->is_hanging(solid_p_index))
912  {
913  nod_pt->unpin(solid_p_index);
914  }
915  }
916  }
917 
918 
920  unsigned nvertex_node() const
921  {
923  }
924 
926  Node* vertex_node_pt(const unsigned& j) const
927  {
929  }
930 
933  unsigned nrecovery_order()
934  {
935  return 2;
936  } // NNODE_1D-1;}
937 
938 
942  Node* interpolating_node_pt(const unsigned& n, const int& value_id)
943 
944  {
945 #ifdef PARANOID
947 #endif
948  // If we are at the value return the solid pressure node
949  if (value_id == 0)
950  {
951  return this->solid_pressure_node_pt(n);
952  }
953  // Otherwise return the nodal values
954  else
955  {
956  return this->node_pt(n);
957  }
958  }
959 
962  double local_one_d_fraction_of_interpolating_node(const unsigned& n1d,
963  const unsigned& i,
964  const int& value_id)
965  {
966 #ifdef PARANOID
968 #endif
969  // If it's the only value, we have the pressure
970  if (value_id == 0)
971  {
972  // The pressure nodes are just located on the boundaries at 0 or 1
973  return double(n1d);
974  }
975  // Otherwise we have the geometric nodes
976  else
977  {
978  return this->local_one_d_fraction_of_node(n1d, i);
979  }
980  }
981 
987  const int& value_id)
988  {
989 #ifdef PARANOID
991 #endif
992 
993  // If we are calculating solid pressure nodes
994  if (value_id == 0)
995  {
996  // Storage for the index of the pressure node
997  unsigned total_index = 0;
998  // The number of nodes along each 1d edge is 2.
999  unsigned NNODE_1D = 2;
1000  // Storage for the index along each boundary
1001  Vector<int> index(DIM);
1002  // Loop over the coordinates
1003  for (unsigned i = 0; i < DIM; i++)
1004  {
1005  // If we are at the lower limit, the index is zero
1006  if (s[i] == -1.0)
1007  {
1008  index[i] = 0;
1009  }
1010  // If we are at the upper limit, the index is the number of nodes
1011  // minus 1
1012  else if (s[i] == 1.0)
1013  {
1014  index[i] = NNODE_1D - 1;
1015  }
1016  // Otherwise, we have to calculate the index in general
1017  else
1018  {
1019  // For uniformly spaced nodes the 0th node number would be
1020  double float_index = 0.5 * (1.0 + s[i]) * (NNODE_1D - 1);
1021  index[i] = int(float_index);
1022  // What is the excess. This should be safe because the
1023  // taking the integer part rounds down
1024  double excess = float_index - index[i];
1025  // If the excess is bigger than our tolerance there is no node,
1026  // return null
1027  if ((excess > FiniteElement::Node_location_tolerance) &&
1028  ((1.0 - excess) > FiniteElement::Node_location_tolerance))
1029  {
1030  return 0;
1031  }
1032  }
1034  total_index +=
1035  index[i] * static_cast<unsigned>(pow(static_cast<float>(NNODE_1D),
1036  static_cast<int>(i)));
1037  }
1038  // If we've got here we have a node, so let's return a pointer to it
1039  return this->solid_pressure_node_pt(total_index);
1040  }
1041  // Otherwise velocity nodes are the same as pressure nodes
1042  else
1043  {
1044  return this->get_node_at_local_coordinate(s);
1045  }
1046  }
1047 
1048 
1051  unsigned ninterpolating_node_1d(const int& value_id)
1052  {
1053 #ifdef PARANOID
1054  RefineableElement::check_value_id(1, value_id);
1055 #endif
1056 
1057  if (value_id == 0)
1058  {
1059  return 2;
1060  }
1061  else
1062  {
1063  return this->nnode_1d();
1064  }
1065  }
1066 
1069  unsigned ninterpolating_node(const int& value_id)
1070  {
1071 #ifdef PARANOID
1072  RefineableElement::check_value_id(1, value_id);
1073 #endif
1074 
1075  if (value_id == 0)
1076  {
1077  return static_cast<unsigned>(pow(2.0, static_cast<int>(DIM)));
1078  }
1079  else
1080  {
1081  return this->nnode();
1082  }
1083  }
1084 
1088  Shape& psi,
1089  const int& value_id) const
1090  {
1091 #ifdef PARANOID
1092  RefineableElement::check_value_id(1, value_id);
1093 #endif
1094 
1095  if (value_id == 0)
1096  {
1097  return this->solid_pshape(s, psi);
1098  }
1099  else
1100  {
1101  return this->shape(s, psi);
1102  }
1103  }
1104 
1105 
1109  {
1110  this->setup_hang_for_value(this->solid_p_nodal_index());
1111  }
1112 
1113  // Return a pointer to the solid node at which pressure dof l2 is stored
1114  Node* solid_pressure_node_pt(const unsigned& l)
1115  {
1116  return this->node_pt(this->Pconv[l]);
1117  }
1118  };
1119 
1120 
1121  //=========================================================================
1124  //=========================================================================
1125  template<>
1127  : public virtual SolidQElement<1, 3>
1128  {
1129  public:
1130  // Make sure that we call the constructor of the SolidQElement
1131  // Only the Intel compiler seems to need this!
1133  };
1134 
1135  //=========================================================================
1138  //=========================================================================
1139  template<>
1142  : public virtual PointElement
1143  {
1144  public:
1145  // Make sure that we call the constructor of the SolidQElement
1146  // Only the Intel compiler seems to need this!
1148  };
1149 
1150  //=========================================================================
1152  //=========================================================================
1153  template<>
1155  : public virtual SolidQElement<2, 3>
1156  {
1157  public:
1158  // Make sure that we call the constructor of the SolidQElement
1159  // Only the Intel compiler seems to need this!
1161  };
1162 
1163  //=========================================================================
1166  //=========================================================================
1167  template<>
1170  : public virtual SolidQElement<1, 3>
1171  {
1172  public:
1173  // Make sure that we call the constructor of the SolidQElement
1174  // Only the Intel compiler seems to need this!
1176  };
1177 
1178 } // namespace oomph
1179 
1180 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_solid_elements.h:1132
FaceGeometry()
Definition: refineable_solid_elements.h:1160
FaceGeometry()
Definition: refineable_solid_elements.h:510
FaceGeometry()
Definition: refineable_solid_elements.h:797
FaceGeometry()
Definition: refineable_solid_elements.h:230
FaceGeometry()
Definition: refineable_solid_elements.h:257
Definition: elements.h:4998
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
static const double Node_location_tolerance
Definition: elements.h:1374
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual Node * get_node_at_local_coordinate(const Vector< double > &s) const
Definition: elements.cc:3882
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
virtual unsigned nnode_1d() const
Definition: elements.h:2218
virtual double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
Definition: elements.h:1858
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
Definition: mesh.h:67
Definition: nodes.h:906
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
static Vector< Vector< int > > Direction_to_vector
Definition: octree.h:353
Definition: oomph_definitions.h:222
bool Unsteady
Flag that switches inertia on/off.
Definition: solid_elements.h:421
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
Definition: solid_elements.h:424
bool is_jacobian_evaluated_by_fd() const
Return the flag indicating whether the jacobian is evaluated by fd.
Definition: solid_elements.h:386
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
Definition: solid_elements.h:140
ConstitutiveLaw *& constitutive_law_pt()
Return the constitutive law pointer.
Definition: solid_elements.h:101
ConstitutiveLaw * Constitutive_law_pt
Pointer to the constitutive law.
Definition: solid_elements.h:415
bool is_inertia_enabled() const
Access function to flag that switches inertia on/off (const version)
Definition: solid_elements.h:164
double *& lambda_sq_pt()
Access function for pointer to timescale ratio (nondim density)
Definition: solid_elements.h:115
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain) const
Return the strain tensor.
Definition: solid_elements.cc:47
IsotropicGrowthFctPt Isotropic_growth_fct_pt
Pointer to isotropic growth function.
Definition: solid_elements.h:409
bool Evaluate_jacobian_by_fd
Use FD to evaluate Jacobian.
Definition: solid_elements.h:430
IsotropicGrowthFctPt & isotropic_growth_fct_pt()
Access function: Pointer to isotropic growth function.
Definition: solid_elements.h:122
double * Lambda_sq_pt
Timescale ratio (non-dim. density)
Definition: solid_elements.h:418
Definition: solid_elements.h:863
double interpolated_solid_p(const Vector< double > &s)
Return the interpolated_solid_pressure.
Definition: solid_elements.h:1011
bool Incompressible
Boolean to determine whether the solid is incompressible or not.
Definition: solid_elements.h:1237
bool is_incompressible() const
Return whether the material is incompressible.
Definition: solid_elements.h:875
Definition: solid_elements.h:440
Definition: elements.h:3439
Definition: solid_elements.h:1523
unsigned npres_solid() const
Return number of pressure values.
Definition: solid_elements.h:1592
void solid_pshape(const Vector< double > &s, Shape &psi) const
Return the pressure shape functions.
int solid_p_nodal_index() const
Set the value at which the solid pressure is stored in the nodes.
Definition: solid_elements.h:1566
static const unsigned Pconv[]
Definition: solid_elements.h:1545
Definition: solid_elements.h:1332
unsigned npres_solid() const
Return number of pressure values.
Definition: solid_elements.h:1390
double solid_p(const unsigned &l)
Return the lth pressure value.
Definition: solid_elements.h:1378
unsigned P_solid_internal_index
Definition: solid_elements.h:1347
Definition: solid_elements.h:608
Definition: refineable_elements.h:97
static void check_value_id(const int &n_continuously_interpolated_values, const int &value_id)
Definition: refineable_elements.cc:53
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Definition: refineable_solid_elements.h:284
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
No values are interpolated in this element (pure solid)
Definition: refineable_solid_elements.h:306
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 0 (pure solid problem)
Definition: refineable_solid_elements.h:373
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_solid_elements.h:321
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
in strain tensor.
Definition: refineable_solid_elements.h:330
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
No values are interpolated in this element (pure solid)
Definition: refineable_solid_elements.h:314
virtual Node * solid_pressure_node_pt(const unsigned &l)
Definition: refineable_solid_elements.h:379
void fill_in_generic_residual_contribution_pvd_with_pressure(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, const unsigned &flag)
Definition: refineable_solid_elements.cc:676
void further_build()
Pass the generic stuff down to the sons.
Definition: refineable_solid_elements.h:386
void get_mass_matrix_diagonal(Vector< double > &mass_diag)
Definition: refineable_solid_elements.cc:551
RefineablePVDEquationsWithPressure()
Constructor:
Definition: refineable_solid_elements.h:287
Class for Refineable PVD equations.
Definition: refineable_solid_elements.h:47
RefineablePVDEquations()
Constructor.
Definition: refineable_solid_elements.h:50
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_solid_elements.h:80
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_solid_elements.h:88
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
No values are interpolated in this element (pure solid)
Definition: refineable_solid_elements.h:73
virtual Node * solid_pressure_node_pt(const unsigned &l)
Definition: refineable_solid_elements.h:138
void further_build()
Further build function, pass the pointers down to the sons.
Definition: refineable_solid_elements.h:144
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
No values are interpolated in this element (pure solid)
Definition: refineable_solid_elements.h:65
void fill_in_generic_contribution_to_residuals_pvd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Call the residuals including hanging node cases.
Definition: refineable_solid_elements.cc:38
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 0 (pure solid problem)
Definition: refineable_solid_elements.h:130
Definition: refineable_brick_element.h:68
Definition: refineable_solid_elements.h:826
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &value_id)
Definition: refineable_solid_elements.h:962
void rebuild_from_sons(Mesh *&mesh_pt)
Empty rebuild from sons, empty.
Definition: refineable_solid_elements.h:853
unsigned required_nvalue(const unsigned &n) const
Definition: refineable_solid_elements.h:841
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &value_id) const
Definition: refineable_solid_elements.h:1087
unsigned ninterpolating_node(const int &value_id)
Definition: refineable_solid_elements.h:1069
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_solid_elements.h:920
RefineableQPVDElementWithContinuousPressure()
Constructor:
Definition: refineable_solid_elements.h:829
void pin_elemental_redundant_nodal_solid_pressures()
Pin the redundant solid pressure.
Definition: refineable_solid_elements.h:893
Node * solid_pressure_node_pt(const unsigned &l)
Definition: refineable_solid_elements.h:1114
void unpin_elemental_solid_pressure_dofs()
Unpin all pressure dofs.
Definition: refineable_solid_elements.h:879
Node * interpolating_node_pt(const unsigned &n, const int &value_id)
Definition: refineable_solid_elements.h:942
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_solid_elements.h:926
unsigned ninterpolating_node_1d(const int &value_id)
Definition: refineable_solid_elements.h:1051
unsigned nrecovery_order()
Definition: refineable_solid_elements.h:933
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
OK, interpolate the solid pressures.
Definition: refineable_solid_elements.h:856
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &value_id)
Definition: refineable_solid_elements.h:986
unsigned ncont_interpolated_values() const
Number of continuously interpolated values (1) solid pressure.
Definition: refineable_solid_elements.h:847
void further_setup_hanging_nodes()
Definition: refineable_solid_elements.h:1108
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
OK get the time-dependent verion.
Definition: refineable_solid_elements.h:867
Definition: refineable_solid_elements.h:436
void further_setup_hanging_nodes()
Definition: refineable_solid_elements.h:486
RefineableQPVDElementWithPressure()
Constructor:
Definition: refineable_solid_elements.h:451
unsigned nrecovery_order()
Definition: refineable_solid_elements.h:479
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_solid_elements.h:466
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_solid_elements.h:472
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 0 (pure solid problem)
Definition: refineable_solid_elements.h:494
void unpin_elemental_solid_pressure_dofs()
Unpin all solid pressure dofs.
Definition: refineable_solid_elements.h:439
Class for refineable QPVDElement elements.
Definition: refineable_solid_elements.h:181
RefineableQPVDElement()
Constructor:
Definition: refineable_solid_elements.h:184
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_solid_elements.h:203
void further_setup_hanging_nodes()
Definition: refineable_solid_elements.h:217
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_solid_elements.h:197
void rebuild_from_sons(Mesh *&mesh_pt)
Empty rebuild from sons, no need to reconstruct anything here.
Definition: refineable_solid_elements.h:194
unsigned nrecovery_order()
Definition: refineable_solid_elements.h:210
Definition: refineable_elements.h:874
virtual void further_build()
Definition: refineable_elements.h:1014
Definition: Qelements.h:2286
Definition: shape.h:76
Definition: Qelements.h:1742
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
#define DIM
Definition: linearised_navier_stokes_elements.h:44
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
@ RDF
Definition: octree.h:54
@ RUB
Definition: octree.h:52
@ LUF
Definition: octree.h:55
@ LDF
Definition: octree.h:53
@ RDB
Definition: octree.h:50
@ LUB
Definition: octree.h:51
@ RUF
Definition: octree.h:56
@ LDB
Definition: octree.h:49
@ SE
Definition: quadtree.h:57
@ NW
Definition: quadtree.h:58
@ NE
Definition: quadtree.h:59
@ SW
Definition: quadtree.h:56
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2