29 #ifndef OOMPH_AXISYM_POROELASITICTY_FACE_ELEMENTS_HEADER
30 #define OOMPH_AXISYM_POROELASITICTY_FACE_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
48 namespace AxisymmetricPoroelasticityTractionElementHelper
58 unsigned n_dim = load.size();
59 for (
unsigned i = 0;
i < n_dim;
i++)
91 template<
class ELEMENT>
122 const unsigned& intpt,
136 const unsigned& intpt,
163 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(element_pt);
165 if (elem_pt->dim() == 3)
174 throw OomphLibError(
"This flux element will not work correctly "
175 "if nodes are hanging\n",
176 OOMPH_CURRENT_FUNCTION,
177 OOMPH_EXCEPTION_LOCATION);
243 const unsigned&
i)
const
257 void output(std::ostream& outfile,
const unsigned& n_plot)
264 unsigned n_node =
nnode();
271 DShape dpsids(n_node, n_dim - 1);
278 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
299 bulk_pt->interpolated_u(s_bulk, disp);
303 bulk_pt->interpolated_du_dt(s_bulk, du_dt);
307 bulk_pt->interpolated_q(s_bulk, q);
310 const double permeability = bulk_pt->permeability();
323 for (
unsigned l = 0; l < n_node; l++)
326 for (
unsigned i = 0;
i < 2;
i++)
329 unsigned u_nodal_index = bulk_pt->u_index_axisym_poroelasticity(
i);
332 interpolated_T1[
i] +=
340 sqrt(1.0 + (interpolated_t1[0] * interpolated_t1[0]) /
341 (interpolated_t1[1] * interpolated_t1[1])) *
346 double J_def = sqrt(1.0 + (interpolated_T1[0] * interpolated_T1[0]) /
347 (interpolated_T1[1] * interpolated_T1[1])) *
362 double lagr_euler_translation_factor =
366 for (
unsigned i = 0;
i < n_dim;
i++)
368 outfile << x[
i] <<
" ";
372 for (
unsigned i = 0;
i < n_dim;
i++)
374 outfile << disp[
i] <<
" ";
378 for (
unsigned i = 0;
i < n_dim;
i++)
387 outfile << permeability * q[0] <<
" "
388 << permeability * q[1] <<
" ";
391 outfile << du_dt[0] <<
" "
395 outfile << du_dt[0] + permeability * q[0] <<
" "
396 << du_dt[1] + permeability * q[1] <<
" ";
399 outfile << unit_normal[0] <<
" "
400 << unit_normal[1] <<
" ";
403 outfile << bulk_pt->interpolated_div_q(s_bulk) <<
" ";
406 outfile << bulk_pt->interpolated_p(s_bulk) <<
" ";
409 outfile << J_undef <<
" ";
412 outfile << J_def <<
" ";
415 outfile << lagr_euler_translation_factor <<
" ";
417 outfile << std::endl;
431 void output(FILE* file_pt,
const unsigned& n_plot)
453 unsigned n_node =
nnode();
459 DShape dpsids(n_node, n_dim - 1);
476 bulk_pt->interpolated_u(s_bulk, disp);
481 for (
unsigned l = 0; l < n_node; l++)
484 for (
unsigned i = 0;
i < 2;
i++)
487 unsigned u_nodal_index = bulk_pt->u_index_axisym_poroelasticity(
i);
490 interpolated_T1[
i] +=
497 double J_undef = sqrt(interpolated_t1[0] * interpolated_t1[0] +
498 interpolated_t1[1] * interpolated_t1[1]) *
503 double J_def = sqrt(interpolated_T1[0] * interpolated_T1[0] +
504 interpolated_T1[1] * interpolated_T1[1]) *
507 double return_val = 1.0;
508 if (J_def != 0.0) return_val = J_undef / J_def;
531 double& seepage_flux_contrib)
537 const double permeability = bulk_el_pt->permeability();
540 const unsigned n_node = this->
nnode();
554 skeleton_flux_contrib = 0.0;
555 seepage_flux_contrib = 0.0;
556 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
570 for (
unsigned l = 0; l < n_node; l++)
573 for (
unsigned i = 0;
i < 2;
i++)
581 double tlength = interpolated_t1[0] * interpolated_t1[0] +
582 interpolated_t1[1] * interpolated_t1[1];
596 bulk_el_pt->interpolated_q(s_bulk, q);
600 bulk_el_pt->interpolated_du_dt(s_bulk, du_dt);
604 x_bulk[0] = bulk_el_pt->interpolated_x(s_bulk, 0);
605 x_bulk[1] = bulk_el_pt->interpolated_x(s_bulk, 1);
609 double tol = 1.0e-10;
612 std::stringstream junk;
613 junk <<
"Gap between bulk and face element coordinate\n"
614 <<
"is suspiciously large: " << error
615 <<
"\nBulk at: " << x_bulk[0] <<
" " << x_bulk[1] <<
"\n"
619 junk.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
625 double dudt_flux = 0.0;
626 for (
unsigned i = 0;
i < 2;
i++)
628 q_flux += permeability * q[
i] * interpolated_normal[
i];
629 dudt_flux += du_dt[
i] * interpolated_normal[
i];
633 seepage_flux_contrib +=
635 skeleton_flux_contrib +=
650 template<
class ELEMENT>
654 unsigned n_dim = this->nodal_dimension();
658 interpolated_x(
s, x);
662 outer_unit_normal(
s, unit_normal);
668 get_traction(time, ipt, x, unit_normal, traction);
676 template<
class ELEMENT>
680 unsigned n_dim = this->nodal_dimension();
684 interpolated_x(
s, x);
688 outer_unit_normal(
s, unit_normal);
694 get_pressure(time, ipt, x, unit_normal, pressure);
702 template<
class ELEMENT>
708 unsigned n_node = nnode();
711 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
715 unsigned n_position_type = this->nnodal_position_type();
716 if (n_position_type != 1)
718 throw OomphLibError(
"Poroelasticity equations are not yet implemented "
719 "for more than one position type",
720 OOMPH_CURRENT_FUNCTION,
721 OOMPH_EXCEPTION_LOCATION);
726 unsigned n_dim = this->nodal_dimension();
730 ELEMENT* bulk_el_pt =
dynamic_cast<ELEMENT*
>(bulk_element_pt());
732 unsigned n_q_basis = bulk_el_pt->nq_basis();
733 unsigned n_q_basis_edge = bulk_el_pt->nq_basis_edge();
742 DShape dpsids(n_node, n_dim - 1);
743 Shape q_basis(n_q_basis, n_dim);
746 unsigned n_intpt = integral_pt()->nweight();
752 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
755 double w = integral_pt()->weight(ipt);
758 dshape_local_at_knot(ipt, psi, dpsids);
761 for (
unsigned i = 0;
i < n_dim - 1;
i++)
763 s_face[
i] = integral_pt()->knot(ipt,
i);
765 s_bulk = local_coordinate_in_bulk(s_face);
768 ELEMENT* bulk_el_pt =
dynamic_cast<ELEMENT*
>(bulk_element_pt());
773 bulk_el_pt->get_q_basis(s_bulk, q_basis);
782 for (
unsigned l = 0; l < n_node; l++)
785 for (
unsigned i = 0;
i < n_dim;
i++)
788 const double x_local = nodal_position(l,
i);
789 interpolated_x[
i] += x_local * psi(l);
792 for (
unsigned j = 0; j < n_dim - 1; j++)
794 interpolated_A(j,
i) += x_local * dpsids(l, j);
801 for (
unsigned i = 0;
i < n_dim - 1;
i++)
803 for (
unsigned j = 0; j < n_dim - 1; j++)
809 for (
unsigned k = 0; k < n_dim; k++)
811 A(
i, j) += interpolated_A(
i, k) * interpolated_A(j, k);
818 outer_unit_normal(ipt, interpolated_normal);
828 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
832 "Wrong dimension in AxisymmetricPoroelasticityTractionElement",
833 OOMPH_CURRENT_FUNCTION,
834 OOMPH_EXCEPTION_LOCATION);
839 double W = w * sqrt(Adet);
843 get_traction(time, ipt, interpolated_x, interpolated_normal, traction);
847 get_pressure(time, ipt, interpolated_x, interpolated_normal, pressure);
850 for (
unsigned l = 0; l < n_node; l++)
853 for (
unsigned i = 0;
i < n_dim;
i++)
855 local_eqn = this->nodal_local_eqn(
856 l, bulk_el_pt->u_index_axisym_poroelasticity(
i));
861 residuals[local_eqn] -=
862 traction[
i] * psi(l) * interpolated_x[0] *
W;
869 for (
unsigned l = 0; l < n_q_basis_edge; l++)
871 local_eqn = this->nodal_local_eqn(1, bulk_el_pt->q_edge_index(l));
877 for (
unsigned i = 0;
i < n_dim;
i++)
880 residuals[local_eqn] += pressure * q_basis(l,
i) *
881 interpolated_normal[
i] * interpolated_x[0] *
900 template<
class POROELASTICITY_BULK_ELEMENT,
class NAVIER_STOKES_BULK_ELEMENT>
903 POROELASTICITY_BULK_ELEMENT>,
920 const double&
q()
const
935 const unsigned& intpt,
941 NAVIER_STOKES_BULK_ELEMENT* ext_el_pt =
942 dynamic_cast<NAVIER_STOKES_BULK_ELEMENT*
>(
955 ext_el_pt->traction(s_ext, interpolated_normal, traction_nst);
963 for (
unsigned i = 0;
i < (n_dim - 1);
i++)
972 x_bulk[0] = ext_el_pt->interpolated_x(s_ext, 0);
973 x_bulk[1] = ext_el_pt->interpolated_x(s_ext, 1);
975 sqrt((x_local[0] - x_bulk[0]) * (x_local[0] - x_bulk[0]) +
976 (x_local[1] - x_bulk[1]) * (x_local[1] - x_bulk[1]));
977 double tol = 1.0e-10;
980 std::stringstream junk;
981 junk <<
"Gap between external and face element coordinate\n"
982 <<
"is suspiciously large:" << error <<
" ( tol = " << tol
984 <<
"\nExternal/bulk at: " << x_bulk[0] <<
" " << x_bulk[1]
986 <<
"Face at: " << x_local[0] <<
" " << x_local[1] <<
"\n";
988 junk.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
995 const double q_value =
q();
1000 traction[0] = q_value * traction_nst[0];
1001 traction[1] = q_value * traction_nst[1];
1007 const unsigned& intpt,
1013 NAVIER_STOKES_BULK_ELEMENT* ext_el_pt =
1014 dynamic_cast<NAVIER_STOKES_BULK_ELEMENT*
>(
1027 ext_el_pt->traction(s_ext, interpolated_normal, traction_nst);
1030 const double q_value =
q();
1035 pressure = -q_value * (interpolated_normal[0] * traction_nst[0] +
1036 interpolated_normal[1] * traction_nst[1]);
1059 void output(std::ostream& outfile,
const unsigned& n_plot)
1077 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1081 for (
unsigned i = 0;
i < n_dim - 1;
i++)
1094 POROELASTICITY_BULK_ELEMENT* bulk_pt =
1099 const double local_permeability = bulk_pt->permeability();
1103 bulk_pt->interpolated_q(s_bulk,
q);
1106 bulk_pt->interpolated_u(s_bulk, disp);
1118 for (
unsigned i = 0;
i < n_dim;
i++)
1120 outfile << x[
i] <<
" ";
1124 for (
unsigned i = 0;
i < n_dim;
i++)
1126 outfile << disp[
i] <<
" ";
1130 for (
unsigned i = 0;
i < n_dim;
i++)
1139 outfile << local_permeability *
q[0] <<
" "
1140 << local_permeability *
q[1] <<
" ";
1144 bulk_pt->interpolated_du_dt(s_bulk, du_dt);
1145 outfile << du_dt[0] <<
" "
1149 outfile << du_dt[0] + local_permeability *
q[0] <<
" "
1150 << du_dt[1] + local_permeability *
q[1] <<
" ";
1153 outfile << unit_normal[0] <<
" "
1154 << unit_normal[1] <<
" ";
1157 outfile << bulk_pt->interpolated_div_q(s_bulk) <<
" ";
1160 outfile << bulk_pt->interpolated_p(s_bulk) <<
" ";
1161 outfile << std::endl;
1186 template<
class POROELASTICITY_BULK_ELEMENT,
class NAVIER_STOKES_BULK_ELEMENT>
1187 double FSILinearisedAxisymPoroelasticTractionElement<
1188 POROELASTICITY_BULK_ELEMENT,
1189 NAVIER_STOKES_BULK_ELEMENT>::Default_Q_Value = 1.0;
A class for elements that allow the imposition of an applied combined traction and pore fluid pressur...
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the traction vector: Pass number of integration point (dummy), Eulerrian coordinate and normal ve...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, double &pressure) pressure_fct_pt()
Reference to the pressure function pointer.
virtual void get_pressure(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, double &pressure)
Get the pressure value: Pass number of integration point (dummy), Eulerrian coordinate and normal vec...
void output(FILE *file_pt)
C_style output function.
void contribution_to_total_porous_flux(double &skeleton_flux_contrib, double &seepage_flux_contrib)
Compute contributions to integrated porous flux over boundary: q_skeleton = \int \partial u_displ / \...
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
AxisymmetricPoroelasticityTractionElement()
Default constructor.
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function. Arguments: Eulerian coordinate; outer unit normal; applied t...
double lagrangian_eulerian_translation_factor(const Vector< double > &s)
Ratio of lengths of line elements (or annular surface areas) in the undeformed and deformed configura...
void output(std::ostream &outfile)
Output function.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void(* Pressure_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, double &result)
Pointer to an imposed pressure function. Arguments: Eulerian coordinate; outer unit normal; applied p...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
void fill_in_contribution_to_residuals_axisym_poroelasticity_face(Vector< double > &residuals)
Helper function that actually calculates the residuals.
void pressure(const double &time, const Vector< double > &s, double &pressure)
Compute pressure value at specified local coordinate Should only be used for post-processing; ignores...
AxisymmetricPoroelasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void traction(const double &time, const Vector< double > &s, Vector< double > &traction)
Compute traction vector at specified local coordinate Should only be used for post-processing; ignore...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Access function to get source element's local coords for specified interaction index at specified int...
void set_ninteraction(const unsigned &n_interaction)
Set the number of interactions in the element This function is usually called in the specific element...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contributions to the jacobian from all external interaction degrees of freedom (geometr...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
double *& q_pt()
Return a pointer the ratio of stress scales used to non-dimensionalise the fluid and poroelastic equa...
void get_pressure(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, double &pressure)
Get the pore fluid pressure from the neighbouring Navier-Stokes bulk element's stress.
FSILinearisedAxisymPoroelasticTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
double * Q_pt
Pointer to the ratio, , of the stress used to non-dimensionalise the fluid stresses to the stress use...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
static double Default_Q_Value
Static default value for the ratio of stress scales used in the fluid and poroelasticity equations (d...
void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the (combined) traction from the neighbouring Navier-Stokes bulk element's stress.
FSILinearisedAxisymPoroelasticTractionElement()
Default constructor.
const double & q() const
Return the ratio of the stress scales used to non-dimensionalise the fluid and poroelasticity equatio...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – overloaded version – ignores n_plot since fsi elements can only evaluate traction a...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
Vector< double > local_coordinate_in_bulk(const Vector< double > &s) const
Return vector of local coordinates in bulk element, given the local coordinates in this FaceElement.
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Return the geometric shape function and its derivative w.r.t. the local coordinates at the ipt-th int...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
unsigned nnode() const
Return the number of nodes.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Function to compute the geometric shape functions and derivatives w.r.t. local coordinates at local c...
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
bool Allow_gap_in_FSI
Public boolean to allow gap between poro-elastic and Navier Stokes element in FSI computations....
void Zero_pressure_fct(const double &time, const Vector< double > &x, const Vector< double > &N, double &load)
Default load function (zero pressure)
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
void output()
Doc the command line arguments.
const double Pi
50 digits from maple
//////////////////////////////////////////////////////////////////// ////////////////////////////////...