28 #ifndef OOMPH_FACE_MESH_PROJECT_HEADER
29 #define OOMPH_FACE_MESH_PROJECT_HEADER
43 template<
class ELEMENT>
57 const unsigned&
i)
const
65 std::ostringstream error_message;
66 error_message <<
"Boundary_id is (still) UINT_MAX -- please set\n"
67 <<
"the actual value with set_boundary_id(...)\n";
93 std::ostringstream error_message;
94 error_message <<
"Boundary_id is (still) UINT_MAX -- please set\n"
95 <<
"the actual value with set_boundary_id(...)\n";
110 unsigned nnod = this->
nnode();
114 for (
unsigned j = 0;
j < nnod;
j++)
117 data_values[
j] = std::make_pair(this->
node_pt(
j), fld);
165 unsigned n_node = this->
nnode();
172 double interpolated_u = 0.0;
175 for (
unsigned l = 0; l < n_node; l++)
177 interpolated_u += this->
nodal_value(t, l, fld) * psi[l];
179 return interpolated_u;
185 return this->
nnode();
221 template<
class GEOMETRIC_ELEMENT>
232 const unsigned& boundary_id,
233 const unsigned& id_of_field_to_be_projected = UINT_MAX)
243 for (
unsigned e = 0;
e < nel;
e++)
250 if (
dynamic_cast<GEOMETRIC_ELEMENT*
>(mesh_pt->
element_pt(
e)) == 0)
252 std::ostringstream error_message;
253 error_message <<
"Element is of wrong type " <<
typeid(el_pt).
name()
254 <<
" doesn't match template parameter!" << std::endl;
261 std::ostringstream error_message;
263 <<
"Internal data will NOT be projected across!\n"
264 <<
"If you want this functionality you'll have to \n"
265 <<
"implement it yourself" << std::endl;
289 unsigned nnod = el_pt->
nnode();
290 for (
unsigned j = 0;
j < nnod;
j++)
295 Node* new_nod_pt = 0;
302 std::ostringstream error_message;
303 error_message <<
"Node isn't on a boundary!" << std::endl;
312 unsigned nval = old_node_pt->
nvalue();
313 unsigned first_index_in_old_node = 0;
317 ->nvalue_assigned_by_face_element(
319 first_index_in_old_node =
321 ->index_of_first_value_assigned_by_face_element(
333 unsigned n_time = old_node_pt->
ntstorage();
334 for (
unsigned t = 0;
t < n_time;
t++)
336 for (
unsigned i = 0;
i < nval;
i++)
339 t,
i, old_node_pt->
value(
t, first_index_in_old_node +
i));
344 unsigned n_dim = old_node_pt->
ndim();
345 for (
unsigned i = 0;
i < n_dim;
i++)
347 new_nod_pt->
x(
i) = old_node_pt->
x(
i);
357 std::ostringstream error_message;
358 error_message <<
"Boundary ID specified as " <<
Boundary_id
359 <<
" but node isn't actually on that boundary!"
408 proj_problem_pt->mesh_pt() = projectable_new_mesh_pt;
411 bool dont_project_positions =
true;
412 proj_problem_pt->project(
this, dont_project_positions);
418 delete proj_problem_pt;
419 delete projectable_new_mesh_pt;
428 for (std::map<Node*, Node*>::iterator it =
New_node_pt.begin();
433 Node* old_node_pt = (*it).first;
436 Node* new_node_pt = (*it).second;
439 unsigned nval = old_node_pt->
nvalue();
440 unsigned first_index_in_old_node = 0;
446 first_index_in_old_node =
448 ->index_of_first_value_assigned_by_face_element(
453 unsigned n_time = old_node_pt->
ntstorage();
454 for (
unsigned t = 0;
t < n_time;
t++)
456 for (
unsigned i = 0;
i < nval;
i++)
459 t, first_index_in_old_node +
i, new_node_pt->
value(
t,
i));
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.)
Definition: face_mesh_project.h:223
BackupMeshForProjection(Mesh *mesh_pt, const unsigned &boundary_id, const unsigned &id_of_field_to_be_projected=UINT_MAX)
Definition: face_mesh_project.h:230
std::map< Node *, Node * > New_node_pt
Map returning new node, labeled by node point in original mesh.
Definition: face_mesh_project.h:468
void copy_onto_original_mesh()
Definition: face_mesh_project.h:426
unsigned Boundary_id
Boundary id.
Definition: face_mesh_project.h:471
unsigned ID_of_field_to_be_projected
Definition: face_mesh_project.h:475
void project_onto_new_mesh(Mesh *new_mesh_pt)
Definition: face_mesh_project.h:394
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
unsigned ntstorage() const
Definition: nodes.cc:879
Definition: elements.h:1313
void set_nodal_dimension(const unsigned &nodal_dim)
Definition: elements.h:1390
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 shape(const Vector< double > &s, Shape &psi) const =0
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 double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
unsigned ninternal_data() const
Return the number of internal data objects.
Definition: elements.h:823
Definition: face_mesh_project.h:46
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: face_mesh_project.h:107
unsigned Boundary_id
Boundary id.
Definition: face_mesh_project.h:199
unsigned nhistory_values_for_coordinate_projection()
Definition: face_mesh_project.h:141
unsigned boundary_id() const
Boundary id.
Definition: face_mesh_project.h:88
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: face_mesh_project.h:160
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Nodal value of boundary coordinate.
Definition: face_mesh_project.h:55
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: face_mesh_project.h:133
GenericLagrangeInterpolatedProjectableElement()
Constructor.
Definition: face_mesh_project.h:49
int local_equation(const unsigned &fld, const unsigned &j)
Definition: face_mesh_project.h:191
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: face_mesh_project.h:183
unsigned nfields_for_projection()
Number of fields to be projected.
Definition: face_mesh_project.h:125
void set_boundary_id(const unsigned &boundary_id)
Boundary id.
Definition: face_mesh_project.h:82
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: face_mesh_project.h:149
Vector< Node * > Node_pt
Vector of pointers to nodes.
Definition: mesh.h:183
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
virtual unsigned ncoordinates_on_boundary(const unsigned &b)
Definition: nodes.cc:2363
virtual bool is_on_boundary() const
Definition: nodes.h:1373
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
virtual void add_to_boundary(const unsigned &b)
Definition: nodes.cc:2336
virtual void set_coordinates_on_boundary(const unsigned &b, const unsigned &k, const Vector< double > &boundary_zeta)
Definition: nodes.cc:2394
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
virtual void get_coordinates_on_boundary(const unsigned &b, const unsigned &k, Vector< double > &boundary_zeta)
Definition: nodes.cc:2379
unsigned nposition_type() const
Definition: nodes.h:1016
double value(const unsigned &i) const
Definition: nodes.cc:2408
Definition: oomph_definitions.h:222
Definition: oomph_definitions.h:267
Definition: projection.h:183
Definition: projection.h:695
unsigned ntstorage() const
Definition: timesteppers.h:601
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
string name
Definition: plotDoE.py:33
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2