29 #ifndef OOMPH_REFINEABLE_ELASTICITY_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_ELASTICITY_ELEMENTS_HEADER
34 #include "../generic/refineable_quad_element.h"
35 #include "../generic/refineable_brick_element.h"
36 #include "../generic/error_estimator.h"
43 template<
unsigned DIM>
62 const unsigned& flag);
92 if (
flux.size() != num_entries)
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
112 for (
unsigned i = 0;
i <
DIM;
i++)
114 flux[icount] = strain(
i,
i);
119 for (
unsigned i = 0;
i <
DIM;
i++)
121 for (
unsigned j =
i + 1;
j <
DIM;
j++)
123 flux[icount] = strain(
i,
j);
177 template<
unsigned DIM,
unsigned NNODE_1D>
223 template<
unsigned NNODE_1D>
236 template<
unsigned NNODE_1D>
250 template<
unsigned NNODE_1D>
263 template<
unsigned NNODE_1D>
279 template<
unsigned DIM>
303 const unsigned& flag);
335 if (
flux.size() != num_entries)
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
355 for (
unsigned i = 0;
i <
DIM;
i++)
357 flux[icount] = strain(
i,
i);
362 for (
unsigned i = 0;
i <
DIM;
i++)
364 for (
unsigned j =
i + 1;
j <
DIM;
j++)
366 flux[icount] = strain(
i,
j);
431 template<
unsigned DIM>
443 for (
unsigned l = 0; l < n_pres; l++)
535 using namespace QuadTreeNames;
538 double centre_solid_p[4] = {0.0, 0.0, 0.0, 0.0};
541 for (
unsigned ison = 0; ison < 4; ison++)
544 centre_solid_p[ison] =
546 quadtree_pt()->son_pt(ison)->object_pt())
551 double p_value = 0.25 * (centre_solid_p[0] + centre_solid_p[1] +
552 centre_solid_p[2] + centre_solid_p[3]);
554 set_solid_p(0, p_value);
564 double slope1 = centre_solid_p[
SE] - centre_solid_p[
SW];
565 double slope2 = centre_solid_p[
NE] - centre_solid_p[
NW];
569 p_value = 0.5 * (slope1 + slope2);
570 set_solid_p(1, p_value);
580 slope1 = centre_solid_p[
NE] - centre_solid_p[
SE];
581 slope2 = centre_solid_p[
NW] - centre_solid_p[
SW];
584 p_value = 0.5 * (slope1 + slope2);
585 set_solid_p(2, p_value);
598 using namespace QuadTreeNames;
601 int son_type = quadtree_pt()->son_type();
617 else if (son_type ==
SE)
623 else if (son_type ==
NE)
630 else if (son_type ==
NW)
642 set_solid_p(0, press);
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));
657 using namespace OcTreeNames;
660 double centre_solid_p[8] = {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0};
663 for (
unsigned ison = 0; ison < 8; ison++)
665 centre_solid_p[ison] = octree_pt()
668 ->internal_data_pt(this->P_solid_internal_index)
678 double av_press = 0.0;
681 for (
unsigned ison = 0; ison < 8; ison++)
683 av_press += centre_solid_p[ison];
687 internal_data_pt(this->P_solid_internal_index)
688 ->set_value(0, 0.125 * av_press);
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];
705 internal_data_pt(this->P_solid_internal_index)
706 ->set_value(1, 0.25 * (slope1 + slope2 + slope3 + slope4));
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];
722 internal_data_pt(this->P_solid_internal_index)
723 ->set_value(2, 0.25 * (slope1 + slope2 + slope3 + slope4));
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];
739 internal_data_pt(this->P_solid_internal_index)
740 ->set_value(3, 0.25 * (slope1 + slope2 + slope3 + slope4));
753 using namespace OcTreeNames;
756 int son_type = octree_pt()->son_type();
760 octree_pt()->father_pt()->object_pt());
765 for (
unsigned i = 0;
i < 3;
i++)
778 set_solid_p(0, press);
781 for (
unsigned i = 1;
i < 4;
i++)
783 set_solid_p(
i, 0.5 * cast_father_element_pt->
solid_p(
i));
821 template<
unsigned DIM>
883 unsigned n_node = this->
nnode();
885 for (
unsigned n = 0;
n < n_node;
n++)
898 unsigned n_node = this->
nnode();
899 for (
unsigned l = 0; l < n_node; l++)
908 for (
unsigned l = 0; l < n_solid_pres; l++)
913 nod_pt->
unpin(solid_p_index);
997 unsigned total_index = 0;
999 unsigned NNODE_1D = 2;
1003 for (
unsigned i = 0;
i <
DIM;
i++)
1012 else if (
s[
i] == 1.0)
1014 index[
i] = NNODE_1D - 1;
1020 double float_index = 0.5 * (1.0 +
s[
i]) * (NNODE_1D - 1);
1021 index[
i] =
int(float_index);
1024 double excess = float_index - index[
i];
1035 index[
i] *
static_cast<unsigned>(
pow(
static_cast<float>(NNODE_1D),
1036 static_cast<int>(
i)));
1077 return static_cast<unsigned>(
pow(2.0,
static_cast<int>(
DIM)));
1081 return this->
nnode();
1089 const int& value_id)
const
1101 return this->
shape(s, psi);
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:1147
FaceGeometry()
Definition: refineable_solid_elements.h:1175
FaceGeometry()
Definition: refineable_solid_elements.h:525
FaceGeometry()
Definition: refineable_solid_elements.h:812
FaceGeometry()
Definition: refineable_solid_elements.h:243
FaceGeometry()
Definition: refineable_solid_elements.h:270
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
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 rebuild_from_sons(Mesh *&mesh_pt)
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: Qelements.h:1742
RealScalar s
Definition: level1_cplx_impl.h:130
#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