29#ifndef OOMPH_AXISYM_FLUID_TRACTION_ELEMENTS_HEADER
30#define OOMPH_AXISYM_FLUID_TRACTION_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
50 namespace AxisymmetricNavierStokesTractionElementHelper
56 const Vector<double>& x,
57 const Vector<double>& N,
58 Vector<double>& load);
70 template<
class ELEMENT>
96 const unsigned&
intpt,
129 this->U_index_axisymmetric_nst_traction.resize(
n_dim + 1);
130 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
132 this->U_index_axisymmetric_nst_traction[
i] =
174 const unsigned&
i)
const
193 return 2 * (n_dim + 1);
200 const unsigned&
nplot)
const
233 for (
unsigned i = 0;
i <
n_dim;
i++)
242 for (
unsigned i = 0;
i <
n_dim;
i++)
270 <<
"Axisymmetric Fluid Traction Navier-Stokes Elements only store "
271 << 2 * (
n_dim + 1) <<
" fields " << std::endl;
302 <<
"Axisymmetric Fluid Traction Navier-Stokes Elements only store "
303 << 2 * (
n_dim + 1) <<
" fields " << std::endl;
346 for (
unsigned i = 0;
i <
n_dim;
i++)
355 for (
unsigned i = 0;
i <
n_dim;
i++)
369 for (
unsigned i = 0;
i <
n_dim;
i++)
375 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
381 for (
unsigned i = 0;
i <
n_dim;
i++)
419 template<
class ELEMENT>
445 template<
class ELEMENT>
459 if (n_position_type != 1)
461 throw OomphLibError(
"AxisymmetricNavierStokes is not yet implemented for "
462 "more than one position type",
473 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
509 for (
unsigned i = 0;
i <
n_dim;
i++)
516 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
525 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
527 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
533 for (
unsigned k = 0;
k <
n_dim;
k++)
552 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
556 "Wrong dimension in AxisymmetricNavierStokesTractionElement",
557 "AxisymmetricNavierStokesTractionElement::fill_in_contribution_to_"
574 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
600 namespace LinearisedFSIAxisymmetricNStNoSlipBCHelper
614 template<
class FLUID_BULK_ELEMENT,
class SOLID_BULK_ELEMENT>
629 const unsigned&
id = 0);
704 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
718 ext_el_pt->interpolated_du_dt_axisymmetric_linear_elasticity(
s_ext,
724 <<
dudt[1] <<
" " <<
dudt[2] <<
" " <<
zeta[0] << std::endl;
798 const unsigned&
flag);
826 template<
class FLUID_BULK_ELEMENT,
class SOLID_BULK_ELEMENT>
829 LinearisedFSIAxisymmetricNStNoSlipBCElementElement(
831 const int& face_index,
859 for (
unsigned i = 0;
i < 3;
i++)
862 dynamic_cast<FLUID_BULK_ELEMENT*
>(
bulk_el_pt)->u_index_axi_nst(
i);
880 template<
class FLUID_BULK_ELEMENT,
class SOLID_BULK_ELEMENT>
883 fill_in_generic_residual_contribution_fsi_no_slip_axisym(
886 const unsigned&
flag)
889 const unsigned n_node = nnode();
895 const unsigned n_intpt = integral_pt()->nweight();
911 for (
unsigned i = 0;
i < (Dim - 1);
i++)
913 s[
i] = integral_pt()->knot(
ipt,
i);
917 double w = integral_pt()->weight(
ipt);
941 bnod_pt->index_of_first_value_assigned_by_face_element(Id);
944 for (
unsigned i = 0;
i < Dim + 1;
i++)
966 for (
unsigned i = 0;
i < Dim + 1;
i++)
972 local_eqn = nodal_local_eqn(
l, U_index_fsi_no_slip_axisym[
i]);
994 bnod_pt->index_of_first_value_assigned_by_face_element(Id) +
1014 l,
bnod_pt->index_of_first_value_assigned_by_face_element(Id) +
i);
1030 nodal_local_eqn(
l2, U_index_fsi_no_slip_axisym[
i]);
A class for elements that allow the imposition of an applied traction in the axisym Navier Stokes eqn...
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 fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
void scalar_value_paraview(std::ofstream &file_out, const unsigned &k, const unsigned &nplot) const
Write values of the k-th scalar field at the plot points. Needs to be implemented for each new specif...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
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...
Vector< unsigned > U_index_axisymmetric_nst_traction
Index at which the i-th velocity component is stored.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
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...
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) const
Get the traction vector: Pass number of integration point (dummy), Eulerian coordinate and normal vec...
AxisymmetricNavierStokesTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void output(FILE *file_pt)
C_style output function.
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
void fill_in_contribution_to_residuals_axisymmetric_nst_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
void output(std::ostream &outfile)
Output function.
A class that contains the information required by Nodes that are located on Mesh boundaries....
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...
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...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
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,...
void add_additional_values(const Vector< unsigned > &nadditional_values, const unsigned &id)
Helper function adding additional values for the unknowns associated with the FaceElement....
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...
double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s....
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point O...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot. Broken virtual; can b...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
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...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
void interpolated_zeta(const Vector< double > &s, Vector< double > &zeta) const
Calculate the interpolated value of zeta, the intrinsic coordinate of the element when viewed as a co...
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned nnode() const
Return the number of nodes.
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 ...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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 shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
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.
A class for elements that allow the imposition of the linearised FSI no slip condition from an adjace...
Vector< unsigned > U_index_fsi_no_slip_axisym
The index at which the unknowns are stored at the nodes.
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
double * St_pt
Pointer to fluid Strouhal number.
void output(FILE *file_pt)
C-style output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: Output at Gauss points; n_plot is ignored.
LinearisedFSIAxisymmetricNStNoSlipBCElementElement(const LinearisedFSIAxisymmetricNStNoSlipBCElementElement &dummy)=delete
Broken copy constructor.
void fill_in_generic_residual_contribution_fsi_no_slip_axisym(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the contri...
double *& st_pt()
Broken assignment operator.
double st() const
Access function for the fluid Strouhal number.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and its Jacobian matrix.
unsigned Dim
The spatial dimension of the problem.
void output(std::ostream &outfile)
Output function.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
An OomphLibError object which should be thrown when an run-time error is encountered....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
double Default_strouhal_number
Default for fluid Strouhal number.
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).