30 #ifndef OOMPH_SPACETIME_FLUID_TRACTION_ELEMENTS_HEADER
31 #define OOMPH_SPACETIME_FLUID_TRACTION_ELEMENTS_HEADER
35 #include <oomph-lib-config.h>
51 template<
class ELEMENT>
52 class NavierStokesSpaceTimeTractionElement
53 :
public virtual FaceGeometry<ELEMENT>,
54 public virtual FaceElement
59 const Vector<double>&
x,
60 const Vector<double>&
n,
61 Vector<double>& result);
69 const bool& called_from_refineable_constructor =
false)
79 if (!called_from_refineable_constructor)
82 if (element_pt->
dim() == 3)
95 std::ostringstream error_message_stream;
99 <<
"This flux element will not work correctly "
100 <<
"if nodes are hanging" << std::endl;
152 residuals, jacobian, 1);
164 void output(std::ostream& outfile,
const unsigned& nplot)
175 const unsigned& flag);
183 const unsigned&
i)
const
226 for (
unsigned i = 0;
i <
Dim - 1;
i++)
236 (*Traction_fct_pt)(time,
x,
n, result);
253 template<
class ELEMENT>
256 Vector<double>& residuals,
258 const unsigned& flag)
261 unsigned n_node = nnode();
270 unsigned n_intpt = integral_pt()->nweight();
276 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
279 double w = integral_pt()->weight(ipt);
283 double J = shape_and_test_at_knot(ipt, psif, testf);
289 double interpolated_t = 0.0;
292 Vector<double> interpolated_x(
Dim - 1, 0.0);
295 Vector<double> interpolated_spacetime_n(
Dim, 0.0);
298 Vector<double> interpolated_n(
Dim - 1, 0.0);
301 Vector<double> traction(
Dim - 1, 0.0);
304 for (
unsigned l = 0; l < n_node; l++)
307 interpolated_t += raw_nodal_position(l,
Dim - 1) * psif(l);
310 for (
unsigned i = 0;
i <
Dim - 1;
i++)
313 interpolated_x[
i] += nodal_position(l,
i) * psif(l);
318 outer_unit_normal(ipt, interpolated_spacetime_n);
323 if (
std::fabs(interpolated_spacetime_n[
Dim - 1]) > 1.0e-15)
326 std::ostringstream error_message_stream;
330 <<
"The component of the outer unit normal in the "
331 <<
"time-direction\nshould be zero but the actual "
332 <<
"value is: " << interpolated_spacetime_n[
Dim - 1] << std::endl;
335 throw OomphLibError(error_message_stream.str(),
342 for (
unsigned i = 0;
i <
Dim - 1;
i++)
345 interpolated_n[
i] = interpolated_spacetime_n[
i];
349 get_traction(interpolated_t, interpolated_x, interpolated_n, traction);
354 for (
unsigned l = 0; l < n_node; l++)
357 for (
unsigned i = 0;
i <
Dim - 1;
i++)
360 local_eqn = u_local_eqn(l,
i);
366 residuals[local_eqn] += traction[
i] * testf[l] *
W;
390 template<
class ELEMENT>
391 class RefineableNavierStokesSpaceTimeTractionElement
392 :
public virtual NavierStokesSpaceTimeTractionElement<ELEMENT>,
393 public virtual NonRefineableElementWithHangingNodes
432 residuals, jacobian, 1);
442 const unsigned& flag);
453 template<
class ELEMENT>
458 const unsigned& flag)
461 unsigned u_nodal_index[this->Dim - 1];
464 for (
unsigned i = 0;
i < this->Dim - 1;
i++)
468 dynamic_cast<ELEMENT*
>(this->bulk_element_pt())->u_index_nst(
i);
472 unsigned n_node = nnode();
481 unsigned n_intpt = integral_pt()->nweight();
487 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
490 double w = integral_pt()->weight(ipt);
494 double J = this->shape_and_test_at_knot(ipt, psif, testf);
500 double interpolated_t = 0.0;
503 Vector<double> interpolated_x(this->Dim - 1, 0.0);
506 Vector<double> interpolated_spacetime_n(this->Dim, 0.0);
509 Vector<double> interpolated_n(this->Dim - 1, 0.0);
512 Vector<double> traction(this->Dim - 1, 0.0);
515 for (
unsigned l = 0; l < n_node; l++)
518 interpolated_t += raw_nodal_position(l, this->Dim - 1) * psif(l);
521 for (
unsigned i = 0;
i < this->Dim - 1;
i++)
524 interpolated_x[
i] += this->nodal_position(l,
i) * psif(l);
529 this->outer_unit_normal(ipt, interpolated_spacetime_n);
534 if (
std::fabs(interpolated_n[this->Dim - 1]) > 1.0e-15)
537 std::ostringstream error_message_stream;
540 error_message_stream <<
"The component of the outer unit normal in the "
541 <<
"time-direction\nshould be zero but the actual "
542 <<
"value is: " << interpolated_n[this->Dim - 1]
546 throw OomphLibError(error_message_stream.str(),
553 for (
unsigned i = 0;
i < this->Dim - 1;
i++)
556 interpolated_n[
i] = interpolated_spacetime_n[
i];
561 interpolated_t, interpolated_x, interpolated_n, traction);
566 unsigned n_master = 1;
569 double hang_weight = 1.0;
572 HangInfo* hang_info_pt = 0;
576 for (
unsigned l = 0; l < n_node; l++)
579 bool is_node_hanging = this->node_pt(l)->is_hanging();
585 hang_info_pt = this->node_pt(l)->hanging_pt();
588 n_master = hang_info_pt->nmaster();
598 for (
unsigned m = 0;
m < n_master;
m++)
601 for (
unsigned i = 0;
i < this->Dim - 1;
i++)
609 local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(
m),
613 hang_weight = hang_info_pt->master_weight(
m);
619 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i]);
629 residuals[local_eqn] += traction[
i] * testf[l] *
W * hang_weight;
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition: BenchSparseUtil.h:23
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
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
Definition: discontinuous_galerkin_equal_order_pressure_spacetime_fluid_traction_elements.h:54
virtual int u_local_eqn(const unsigned &n, const unsigned &i)
Definition: spacetime_fluid_traction_elements.h:194
unsigned Dim
The highest dimension of the problem.
Definition: discontinuous_galerkin_equal_order_pressure_spacetime_fluid_traction_elements.h:300
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: spacetime_fluid_traction_elements.h:137
NavierStokesSpaceTimeTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Definition: spacetime_fluid_traction_elements.h:66
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
Definition: spacetime_fluid_traction_elements.h:164
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: spacetime_fluid_traction_elements.h:201
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
Definition: spacetime_fluid_traction_elements.h:147
void output(std::ostream &outfile)
Overload the output function.
Definition: spacetime_fluid_traction_elements.h:157
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: spacetime_fluid_traction_elements.h:181
~NavierStokesSpaceTimeTractionElement()
Destructor should not delete anything.
Definition: spacetime_fluid_traction_elements.h:122
void fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: discontinuous_galerkin_equal_order_pressure_spacetime_fluid_traction_elements.h:313
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: discontinuous_galerkin_equal_order_pressure_spacetime_fluid_traction_elements.h:57
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: spacetime_fluid_traction_elements.h:217
void(*&)(const double &t, const Vector< double > &x, const Vector< double > &n, Vector< double > &result) traction_fct_pt()
Access function for the imposed traction pointer.
Definition: discontinuous_galerkin_equal_order_pressure_spacetime_fluid_traction_elements.h:125
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
Definition: oomph_definitions.h:222
Definition: refineable_elements.h:97
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
Definition: spacetime_fluid_traction_elements.h:427
RefineableNavierStokesSpaceTimeTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
Definition: spacetime_fluid_traction_elements.h:397
void refineable_fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: discontinuous_galerkin_equal_order_pressure_spacetime_fluid_traction_elements.h:513
~RefineableNavierStokesSpaceTimeTractionElement()
Destructor should not delete anything.
Definition: spacetime_fluid_traction_elements.h:407
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: spacetime_fluid_traction_elements.h:419
unsigned ncont_interpolated_values() const
Definition: spacetime_fluid_traction_elements.h:412
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
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
@ 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