30 #ifndef OOMPH_POLAR_FLUID_TRACTION_ELEMENTS_HEADER
31 #define OOMPH_POLAR_FLUID_TRACTION_ELEMENTS_HEADER
40 #include "../generic/Qelements.h"
51 template<
class ELEMENT>
82 unsigned n_node =
nnode();
86 for (
unsigned i = 0;
i < n_node;
i++)
104 for (
unsigned i = 0;
i <
Dim;
i++)
112 (*Traction_fct_pt)(time,
x, result);
203 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(element_pt);
205 if (elem_pt->dim() == 3)
214 throw OomphLibError(
"This flux element will not work correctly "
215 "if nodes are hanging\n",
287 void output(std::ostream& outfile,
const unsigned& nplot)
293 double u(
const unsigned& l,
const unsigned&
i)
299 double x(
const unsigned& l,
const unsigned&
i)
315 template<
class ELEMENT>
323 unsigned n_node = nnode();
326 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
329 Shape psif(n_node), testf(n_node);
332 unsigned n_intpt = integral_pt()->nweight();
346 const double eta = get_eta();
349 int local_eqn = 0, local_unknown = 0, pext_local_eqn = 0,
350 pext_local_unknown = 0;
365 pext_local_eqn = external_local_eqn(External_data_number_of_Pext, 0);
372 pext_local_unknown = pext_local_eqn;
377 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
380 double w = integral_pt()->weight(ipt);
384 double J = shape_and_test_at_knot(ipt, psif, testf);
394 for (
unsigned i = 0;
i <
Dim;
i++)
396 interpolated_x[
i] = 0.0;
397 interpolated_u[
i] = 0.0;
402 for (
unsigned l = 0; l < n_node; l++)
405 for (
unsigned i = 0;
i <
Dim;
i++)
408 interpolated_u[
i] += this->nodal_value(l,
i) * psif[l];
409 interpolated_x[
i] += this->nodal_position(l,
i) * psif[l];
415 get_traction(time, interpolated_x, traction);
420 for (
unsigned l = 0; l < n_node; l++)
425 local_eqn = u_local_eqn(l,
i);
431 (interpolated_u[
i] / interpolated_x[0]) *
432 testf[l] * interpolated_x[0] *
Alpha *
W;
438 residuals[local_eqn] +=
439 pext * testf[l] * interpolated_x[0] *
Alpha *
W;
447 for (
unsigned l2 = 0; l2 < n_node; l2++)
453 local_unknown = u_local_eqn(l2, i2);
454 if (local_unknown >= 0)
457 jacobian(local_eqn, local_unknown) -=
459 testf[l] * interpolated_x[0] *
Alpha *
W;
469 if (pext_local_unknown >= 0)
472 jacobian(local_eqn, pext_local_unknown) +=
473 testf[l] * interpolated_x[0] *
Alpha *
W;
492 if (pext_local_eqn >= 0)
495 residuals[pext_local_eqn] +=
496 interpolated_u[0] * interpolated_x[0] *
Alpha *
W;
509 for (
unsigned l2 = 0; l2 < n_node; l2++)
515 local_unknown = u_local_eqn(l2, i2);
516 if (local_unknown >= 0)
519 jacobian(pext_local_eqn, local_unknown) +=
520 psif[l2] * interpolated_x[0] *
Alpha *
W;
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
double value(const unsigned &i) const
Definition: nodes.h:293
Definition: elements.h:4338
int & face_index()
Definition: elements.h:4626
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
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
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 nnode() const
Return the number of nodes.
Definition: elements.h:2210
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220
bool has_hanging_nodes() const
Definition: elements.h:2470
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:307
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
Definition: oomph_definitions.h:222
Definition: polar_fluid_traction_elements.h:54
void output(std::ostream &outfile)
Overload the output function.
Definition: polar_fluid_traction_elements.h:281
const int boundary() const
Boundary.
Definition: polar_fluid_traction_elements.h:167
void get_traction(double time, const Vector< double > &x, Vector< double > &result)
Function to calculate the traction applied to the fluid.
Definition: polar_fluid_traction_elements.h:96
int Boundary
Definition: polar_fluid_traction_elements.h:137
double * Alpha_pt
Pointer to the angle alpha.
Definition: polar_fluid_traction_elements.h:124
~PolarNavierStokesTractionElement()
Destructor should not delete anything.
Definition: polar_fluid_traction_elements.h:238
PolarNavierStokesTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: polar_fluid_traction_elements.h:192
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
Definition: polar_fluid_traction_elements.h:260
void(*&)(const double &t, const Vector< double > &x, Vector< double > &result) traction_fct_pt()
Definition: polar_fluid_traction_elements.h:241
void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: polar_fluid_traction_elements.h:317
unsigned External_data_number_of_Pext
Definition: polar_fluid_traction_elements.h:131
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to an imposed traction function.
Definition: polar_fluid_traction_elements.h:57
unsigned Dim
The highest dimension of the problem.
Definition: polar_fluid_traction_elements.h:62
Data * Pext_data_pt
Pointer to the Data item that stores the external pressure.
Definition: polar_fluid_traction_elements.h:127
const double get_eta() const
Eta.
Definition: polar_fluid_traction_elements.h:179
double x(const unsigned &l, const unsigned &i)
local position
Definition: polar_fluid_traction_elements.h:299
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: polar_fluid_traction_elements.h:77
double u(const unsigned &l, const unsigned &i)
local velocities
Definition: polar_fluid_traction_elements.h:293
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: polar_fluid_traction_elements.h:270
void set_external_pressure_data(Data *pext_data_pt)
Function for setting up external pressure.
Definition: polar_fluid_traction_elements.h:156
void set_eta(double eta)
Function to set Eta.
Definition: polar_fluid_traction_elements.h:185
void set_boundary(int bound)
Function to set boundary.
Definition: polar_fluid_traction_elements.h:173
double *& alpha_pt()
Pointer to Alpha.
Definition: polar_fluid_traction_elements.h:150
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: polar_fluid_traction_elements.h:249
virtual int u_local_eqn(const unsigned &n, const unsigned &i)
Definition: polar_fluid_traction_elements.h:70
double Eta
Definition: polar_fluid_traction_elements.h:140
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
Definition: polar_fluid_traction_elements.h:287
const double & alpha() const
Alpha.
Definition: polar_fluid_traction_elements.h:144
Definition: refineable_elements.h:97
RealScalar alpha
Definition: level1_cplx_impl.h:151
double multiplier(const Vector< double > &xi)
Definition: disk_oscillation.cc:63
Data * Pext_data_pt
Pointer to pressure load.
Definition: steady_ring.cc:60
static const unsigned Dim
Problem dimension.
Definition: two_d_tilted_square.cc:62
double Alpha
Parameter for steepness of step.
Definition: two_d_adv_diff_adapt.cc:53
double eta
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:45
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
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