27 #ifndef OOMPH_NAVIER_STOKES_FLUX_CONTROL_ELEMENTS
28 #define OOMPH_NAVIER_STOKES_FLUX_CONTROL_ELEMENTS
32 #include <oomph-lib-config.h>
36 #include "../generic/nodes.h"
37 #include "../navier_stokes/navier_stokes_surface_power_elements.h"
111 double* prescribed_flux_value_pt)
128 for (
unsigned e = 0;
e < n_el;
e++)
205 std::ostringstream error_message;
206 error_message <<
"Dof_number_for_unknown hasn't been set yet!\n"
207 <<
"Please do so using the dof_number_for_unknown()\n"
208 <<
"access function\n";
245 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
250 std::ostringstream error_message;
251 error_message <<
"Dof_number_for_unknown hasn't been set yet!\n"
252 <<
"Please do so using the dof_number_for_unknown()\n"
253 <<
"access function\n";
261 std::pair<unsigned, unsigned> dof_lookup;
267 dof_lookup_list.push_front(dof_lookup);
277 double volume_flux = 0.0;
281 for (
unsigned e = 0;
e < n_el;
e++)
292 if (flux_control_el_pt == 0)
295 "NavierStokesFluxControlElements",
347 template<
class ELEMENT>
357 const bool& called_from_refineable_constructor =
false)
363 if (!called_from_refineable_constructor)
365 ELEMENT* elem_pt =
new ELEMENT;
367 if (elem_pt->dim() == 3)
373 std::ostringstream error_message;
375 <<
"This element does not work properly with refineable bulk \n"
376 <<
"elements in 3D. Please use the refineable version\n"
412 residuals, jacobian, 1);
439 unsigned n_node = this->
nnode();
443 for (
unsigned i = 0;
i < n_node;
i++)
460 unsigned n_node = this->
nnode();
463 Shape psif(n_node), testf(n_node);
475 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
493 for (
unsigned i = 0;
i <
Dim;
i++)
495 traction[
i] = -pressure * unit_normal[
i];
499 for (
unsigned l = 0; l < n_node; l++)
502 for (
unsigned i = 0;
i <
Dim;
i++)
510 residuals[local_eqn] += traction[
i] * testf[l] *
W;
524 if (local_unknown >= 0)
527 double jac_contribution = -unit_normal[
i] * testf[l] *
W;
528 jacobian(local_eqn, local_unknown) += jac_contribution;
531 jacobian(local_unknown, local_eqn) += jac_contribution;
562 template<
class ELEMENT>
606 residuals, jacobian, 1);
618 unsigned u_nodal_index[this->
Dim];
619 for (
unsigned i = 0;
i < this->
Dim;
i++)
629 unsigned n_node = this->
nnode();
632 Shape psif(n_node), testf(n_node);
645 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
663 for (
unsigned i = 0;
i < this->
Dim;
i++)
665 traction[
i] = -pressure * unit_normal[
i];
671 unsigned n_master = 1;
672 double hang_weight = 1.0;
676 for (
unsigned l = 0; l < n_node; l++)
687 n_master = hang_info_pt->
nmaster();
696 for (
unsigned m = 0;
m < n_master;
m++)
699 for (
unsigned i = 0;
i < this->
Dim;
i++)
725 residuals[local_eqn] +=
726 traction[
i] * testf[l] *
W * hang_weight;
740 if (local_unknown >= 0)
743 double jac_contribution =
744 -unit_normal[
i] * testf[l] *
W * hang_weight;
745 jacobian(local_eqn, local_unknown) += jac_contribution;
748 jacobian(local_unknown, local_eqn) += jac_contribution;
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
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
int & face_index()
Definition: elements.h:4626
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
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:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
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
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220
Definition: elements.h:73
unsigned long eqn_number(const unsigned &ieqn_local) const
Definition: elements.h:704
Data *& external_data_pt(const unsigned &i)
Return a pointer to i-th external data object.
Definition: elements.h:659
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
int external_local_eqn(const unsigned &i, const unsigned &j)
Definition: elements.h:311
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:62
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
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.
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
Definition: navier_stokes_flux_control_elements.h:351
void fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: flux_control_elements_bk.h:324
double get_volume_flux()
Function to get the integral of the volume flux.
Definition: navier_stokes_flux_control_elements.h:416
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: navier_stokes_flux_control_elements.h:434
~NavierStokesFluxControlElement()
Destructor should not delete anything.
Definition: navier_stokes_flux_control_elements.h:393
void fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: navier_stokes_flux_control_elements.h:456
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: navier_stokes_flux_control_elements.h:407
unsigned Dim
The highest dimension of the problem.
Definition: flux_control_elements_bk.h:413
NavierStokesFluxControlElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor, which takes a "bulk" element and face index.
Definition: navier_stokes_flux_control_elements.h:354
virtual int u_local_eqn(const unsigned &n, const unsigned &i)
Definition: navier_stokes_flux_control_elements.h:427
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: navier_stokes_flux_control_elements.h:396
Definition: navier_stokes_surface_power_elements.h:52
double get_volume_flux()
Get integral of volume flux.
Definition: navier_stokes_surface_power_elements.h:704
Definition: flux_control_elements_bk.h:54
double * Prescribed_flux_value_pt
Pointer to the value that stores the prescribed flux.
Definition: navier_stokes_flux_control_elements.h:319
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: navier_stokes_flux_control_elements.h:244
unsigned ndof_types() const
Definition: navier_stokes_flux_control_elements.h:200
unsigned & dof_number_for_unknown()
Definition: navier_stokes_flux_control_elements.h:229
void fill_in_generic_residual_contribution_flux_control(Vector< double > &residuals)
Definition: navier_stokes_flux_control_elements.h:273
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: navier_stokes_flux_control_elements.h:188
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Definition: navier_stokes_flux_control_elements.h:178
unsigned Dim
spatial dim of NS system
Definition: navier_stokes_flux_control_elements.h:329
unsigned Dof_number_for_unknown
Definition: navier_stokes_flux_control_elements.h:326
Data * Pressure_data_pt
Definition: flux_control_elements_bk.h:182
Data * pressure_data_pt() const
Definition: navier_stokes_flux_control_elements.h:170
NetFluxControlElement(Mesh *flux_control_mesh_pt, double *prescribed_flux_value_pt)
Definition: navier_stokes_flux_control_elements.h:110
~NetFluxControlElement()
Empty Destructor - Data gets deleted automatically.
Definition: navier_stokes_flux_control_elements.h:148
unsigned dim() const
Broken assignment operator.
Definition: navier_stokes_flux_control_elements.h:162
Mesh * Flux_control_mesh_pt
Definition: flux_control_elements_bk.h:186
NetFluxControlElement(const NetFluxControlElement &dummy)=delete
Broken copy constructor.
HangInfo *const & hanging_pt() const
Definition: nodes.h:1228
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
Definition: refineable_elements.h:798
Definition: oomph_definitions.h:222
Definition: refineable_elements.h:97
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Definition: refineable_elements.h:278
Definition: navier_stokes_flux_control_elements.h:566
unsigned ncont_interpolated_values() const
Definition: navier_stokes_flux_control_elements.h:583
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: navier_stokes_flux_control_elements.h:601
void refineable_fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: navier_stokes_flux_control_elements.h:614
~RefineableNavierStokesFluxControlElement()
Destructor should not delete anything.
Definition: navier_stokes_flux_control_elements.h:578
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: navier_stokes_flux_control_elements.h:590
RefineableNavierStokesFluxControlElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
Definition: navier_stokes_flux_control_elements.h:569
Definition: navier_stokes_flux_control_elements.h:55
virtual double get_volume_flux()=0
Pure virtual function to calculate integral of the volume flux.
TemplateFreeNavierStokesFluxControlElementBase()
Empty constructor.
Definition: navier_stokes_flux_control_elements.h:58
virtual ~TemplateFreeNavierStokesFluxControlElementBase()
Empty virtual destructor.
Definition: navier_stokes_flux_control_elements.h:61
void add_pressure_data(Data *pressure_data_pt)
Definition: navier_stokes_flux_control_elements.h:68
unsigned & pressure_data_id()
Definition: navier_stokes_flux_control_elements.h:76
unsigned Pressure_data_id
Definition: navier_stokes_flux_control_elements.h:84
int * m
Definition: level2_cplx_impl.h:294
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86