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";
114 for (
unsigned j = 0;
j <
nnod;
j++)
172 double interpolated_u = 0.0;
179 return interpolated_u;
185 return this->
nnode();
221 template<
class GEOMETRIC_ELEMENT>
232 const unsigned& boundary_id,
243 for (
unsigned e = 0;
e <
nel;
e++)
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;
290 for (
unsigned j = 0;
j <
nnod;
j++)
302 std::ostringstream error_message;
303 error_message <<
"Node isn't on a boundary!" << std::endl;
317 ->nvalue_assigned_by_face_element(
321 ->index_of_first_value_assigned_by_face_element(
336 for (
unsigned i = 0;
i <
nval;
i++)
345 for (
unsigned i = 0;
i <
n_dim;
i++)
357 std::ostringstream error_message;
358 error_message <<
"Boundary ID specified as " <<
Boundary_id
359 <<
" but node isn't actually on that boundary!"
428 for (std::map<Node*, Node*>::iterator
it =
New_node_pt.begin();
448 ->index_of_first_value_assigned_by_face_element(
456 for (
unsigned i = 0;
i <
nval;
i++)
Class that makes backup (via a deep copy) of a mesh, keeping alive enough information to allow the so...
BackupMeshForProjection(Mesh *mesh_pt, const unsigned &boundary_id, const unsigned &id_of_field_to_be_projected=UINT_MAX)
Constructor: Pass existing mesh and the boundary ID (need to find the boundary coordinates that are u...
std::map< Node *, Node * > New_node_pt
Map returning new node, labeled by node point in original mesh.
void copy_onto_original_mesh()
Copy nodal values back onto original mesh from which this mesh was built. This obviously only makes s...
unsigned Boundary_id
Boundary id.
unsigned ID_of_field_to_be_projected
ID of field to be projected (UINT_MAX if there isn't one, in which case we're doing all of them.
void project_onto_new_mesh(Mesh *new_mesh_pt)
Project the solution that was present in the original mesh and from which this backup mesh was create...
A class that contains the information required by Nodes that are located on Mesh boundaries....
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
A general Finite Element class.
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
void set_nodal_dimension(const unsigned &nodal_dim)
Set the dimension of the nodes in the element. This will typically only be required when constructing...
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
unsigned ninternal_data() const
Return the number of internal data objects.
Class that makes the finite element specified as template argument projectable – on the assumption th...
unsigned Boundary_id
Boundary id.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!). Extract from first node bu...
unsigned boundary_id() const
Boundary id.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Nodal value of boundary coordinate.
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field (includes current value!). Extract from first ...
GenericLagrangeInterpolatedProjectableElement()
Constructor.
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld. Assumed to be the local nodal equation.
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
unsigned nfields_for_projection()
Number of fields to be projected.
void set_boundary_id(const unsigned &boundary_id)
Boundary id.
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
unsigned ndim() const
Access function to # of Eulerian coordinates.
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
Vector< Node * > Node_pt
Vector of pointers to nodes.
FiniteElement * finite_element_pt(const unsigned &e) const
Upcast (downcast?) to FiniteElement (needed to access FiniteElement member functions).
unsigned long nnode() const
Return number of nodes in the mesh.
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
unsigned long nelement() const
Return number of elements in the mesh.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
virtual void get_coordinates_on_boundary(const unsigned &b, const unsigned &k, Vector< double > &boundary_zeta)
Return the vector of the k-th generalised boundary coordinates on mesh boundary b....
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Projection problem. This is created during the adaptation of unstructured meshes and it is assumed th...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).