30 #ifndef OOMPH_FLUID_TRACTION_ELEMENTS_HEADER
31 #define OOMPH_FLUID_TRACTION_ELEMENTS_HEADER
35 #include <oomph-lib-config.h>
40 #include "../generic/Qelements.h"
41 #include "../generic/Telements.h"
52 template<
class ELEMENT>
69 const unsigned&
i)
const
91 unsigned n_node =
nnode();
95 for (
unsigned i = 0;
i < n_node;
i++)
114 for (
unsigned i = 0;
i <
Dim;
i++)
122 (*Traction_fct_pt)(time,
x,
n, result);
142 const bool& called_from_refineable_constructor =
false)
152 if (!called_from_refineable_constructor)
155 if (element_pt->
dim() == 3)
164 throw OomphLibError(
"This flux element will not work correctly "
165 "if nodes are hanging\n",
209 residuals, jacobian, 1);
219 void output(std::ostream& outfile,
const unsigned& nplot)
235 template<
class ELEMENT>
241 unsigned n_node = nnode();
244 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
247 Shape psif(n_node), testf(n_node);
250 unsigned n_intpt = integral_pt()->nweight();
256 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
259 double w = integral_pt()->weight(ipt);
263 double J = shape_and_test_at_knot(ipt, psif, testf);
272 for (
unsigned i = 0;
i <
Dim;
i++)
274 interpolated_x[
i] = 0.0;
278 for (
unsigned l = 0; l < n_node; l++)
281 for (
unsigned i = 0;
i <
Dim;
i++)
283 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
289 outer_unit_normal(ipt, interpolated_n);
293 get_traction(time, interpolated_x, interpolated_n, traction);
298 for (
unsigned l = 0; l < n_node; l++)
301 for (
unsigned i = 0;
i <
Dim;
i++)
303 local_eqn = u_local_eqn(l,
i);
308 residuals[local_eqn] += traction[
i] * testf[l] *
W;
334 template<
class ELEMENT>
374 residuals, jacobian, 1);
396 template<
class ELEMENT>
402 unsigned u_nodal_index[this->
Dim];
403 for (
unsigned i = 0;
i < this->
Dim;
i++)
406 dynamic_cast<ELEMENT*
>(this->bulk_element_pt())->u_index_nst(
i);
410 unsigned n_node = nnode();
413 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
416 Shape psif(n_node), testf(n_node);
419 unsigned n_intpt = integral_pt()->nweight();
425 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
428 double w = integral_pt()->weight(ipt);
432 double J = this->shape_and_test_at_knot(ipt, psif, testf);
441 for (
unsigned i = 0;
i < this->
Dim;
i++)
443 interpolated_x[
i] = 0.0;
447 for (
unsigned l = 0; l < n_node; l++)
450 for (
unsigned i = 0;
i < this->
Dim;
i++)
452 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
458 this->outer_unit_normal(ipt, interpolated_n);
462 this->get_traction(time, interpolated_x, interpolated_n, traction);
467 unsigned n_master = 1;
468 double hang_weight = 1.0;
475 for (
unsigned l = 0; l < n_node; l++)
478 bool is_node_hanging = this->node_pt(l)->is_hanging();
483 hang_info_pt = this->node_pt(l)->hanging_pt();
486 n_master = hang_info_pt->
nmaster();
495 for (
unsigned m = 0;
m < n_master;
m++)
498 for (
unsigned i = 0;
i < this->
Dim;
i++)
514 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i]);
536 residuals[local_eqn] += traction[
i] * testf[l] *
W * hang_weight;
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Definition: elements.h:4338
int & face_index()
Definition: elements.h:4626
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
double J_eulerian_at_knot(const unsigned &ipt) const
Definition: elements.cc:5328
Definition: elements.h:4998
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220
bool has_hanging_nodes() const
Definition: elements.h:2470
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
Definition: fluid_traction_elements.h:55
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: fluid_traction_elements.h:195
void output(std::ostream &outfile)
Overload the output function.
Definition: fluid_traction_elements.h:213
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function.
Definition: fluid_traction_elements.h:58
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: fluid_traction_elements.h:86
virtual int u_local_eqn(const unsigned &n, const unsigned &i)
Definition: fluid_traction_elements.h:79
void(*&)(const double &t, const Vector< double > &x, const Vector< double > &n, Vector< double > &result) traction_fct_pt()
Definition: fluid_traction_elements.h:186
void fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: fluid_traction_elements.h:237
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
Definition: fluid_traction_elements.h:204
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: fluid_traction_elements.h:67
unsigned Dim
The highest dimension of the problem.
Definition: fluid_traction_elements.h:134
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
Definition: fluid_traction_elements.h:219
~NavierStokesTractionElement()
Destructor should not delete anything.
Definition: fluid_traction_elements.h:183
void get_traction(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Function to calculate the traction applied to the fluid.
Definition: fluid_traction_elements.h:105
NavierStokesTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Definition: fluid_traction_elements.h:139
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
Definition: refineable_elements.h:798
Definition: oomph_definitions.h:222
Definition: refineable_elements.h:97
Definition: fluid_traction_elements.h:338
void refineable_fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: fluid_traction_elements.h:398
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
Definition: fluid_traction_elements.h:369
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: fluid_traction_elements.h:361
unsigned ncont_interpolated_values() const
Definition: fluid_traction_elements.h:354
RefineableNavierStokesTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
Definition: fluid_traction_elements.h:341
~RefineableNavierStokesTractionElement()
Destructor should not delete anything.
Definition: fluid_traction_elements.h:349
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
static const unsigned Dim
Problem dimension.
Definition: two_d_tilted_square.cc:62
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86