Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
oomph::Node Class Reference

Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a given dimension. More...

#include <nodes.h>

+ Inheritance diagram for oomph::Node:

Public Types

typedef void(* AuxNodeUpdateFctPt) (Node *)
 Function pointer to auxiliary node update function. More...
 

Public Member Functions

 Node ()
 Default constructor. More...
 
 Node (const unsigned &n_dim, const unsigned &n_position_type, const unsigned &initial_n_value, const bool &allocate_x_position=true)
 Steady constructor, for a Node of spatial dimension n_dim. Allocates storage for initial_n_value values. NPosition_type is the number of coordinate types needed in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements, etc). More...
 
 Node (TimeStepper *const &time_stepper_pt, const unsigned &n_dim, const unsigned &n_position_type, const unsigned &initial_n_value, const bool &allocate_x_position=true)
 Unsteady constructor for a node of spatial dimension n_dim. Allocates storage for initial_n_value values with history values as required by the timestepper. n_position_type: # of coordinate types needed in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements). More...
 
virtual ~Node ()
 Destructor: Clean up the memory allocated for nodal position. More...
 
 Node (const Node &node)=delete
 Broken copy constructor. More...
 
void operator= (const Node &)=delete
 Broken assignment operator. More...
 
unsigned nposition_type () const
 Number of coordinate types needed in the mapping between local and global coordinates. More...
 
TimeStepper *& position_time_stepper_pt ()
 Return a pointer to the position timestepper. More...
 
TimeStepper *const & position_time_stepper_pt () const
 Return a pointer to the position timestepper (const version). More...
 
virtual void set_position_time_stepper (TimeStepper *const &position_time_stepper_pt, const bool &preserve_existing_data)
 Set a new position timestepper be resizing the appropriate storage. More...
 
virtual bool does_pointer_correspond_to_position_data (double *const &parameter_pt)
 Check whether the pointer parameter_pt addresses position data values. It never does for a standard node, because the positions are not data. More...
 
virtual void assign_eqn_numbers (unsigned long &global_ndof, Vector< double * > &dof_pt)
 Assign global equation numbers; increment global number of unknowns, global_ndof; and add any new dofs to the dof_pt. More...
 
unsigned ndim () const
 Return (Eulerian) spatial dimension of the node. More...
 
double & x (const unsigned &i)
 Return the i-th nodal coordinate. More...
 
const double & x (const unsigned &i) const
 Return the i-th nodal coordinate (const version). More...
 
double & x (const unsigned &t, const unsigned &i)
 Return the position x(i) at previous timestep t (t=0: present; t>0 previous timestep). More...
 
const double & x (const unsigned &t, const unsigned &i) const
 Return the position x(i) at previous timestep t (t=0: present; t>0 previous timestep) (const version) More...
 
double dx_dt (const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt. More...
 
double dx_dt (const unsigned &j, const unsigned &i) const
 Return the i-th component of j-th derivative of nodal position: d^jx/dt^j. More...
 
virtual Nodecopied_node_pt () const
 Return pointer to copied node (null if the current node is not a copy – always the case here; it's overloaded for boundary nodes) More...
 
virtual bool position_is_a_copy () const
 Return whether any position coordinate has been copied (always false) More...
 
virtual bool position_is_a_copy (const unsigned &i) const
 Return whether the position coordinate i has been copied (always false) More...
 
double & x_gen (const unsigned &k, const unsigned &i)
 Reference to the generalised position x(k,i). ‘Type’: k; Coordinate direction: i. More...
 
const double & x_gen (const unsigned &k, const unsigned &i) const
 Reference to the generalised position x(k,i). ‘Type’: k; Coordinate direction: i (const version). More...
 
double & x_gen (const unsigned &t, const unsigned &k, const unsigned &i)
 Reference to the generalised position x(k,i) at the previous timestep [t=0: present]. ‘Type’: k; Coordinate direction: i. More...
 
const double & x_gen (const unsigned &t, const unsigned &k, const unsigned &i) const
 Reference to the generalised position x(k,i) at the previous timestep [t=0: present]. ‘Type’: k; Coordinate direction: i. (const version) More...
 
double dx_gen_dt (const unsigned &k, const unsigned &i) const
 i-th component of time derivative (velocity) of the generalised position, dx(k,i)/dt. ‘Type’: k; Coordinate direction: i. More...
 
double dx_gen_dt (const unsigned &j, const unsigned &k, const unsigned &i) const
 i-th component of j-th time derivative (velocity) of the generalised position, d^jx(k,i)/dt^j. ‘Type’: k; Coordinate direction: i. More...
 
double * x_pt (const unsigned &t, const unsigned &i)
 Direct access to the i-th coordinate at time level t (t=0: present; t>0: previous) More...
 
void copy (Node *orig_node_pt)
 Copy all nodal data from specified Node object. More...
 
virtual void dump (std::ostream &dump_file) const
 Dump nodal position and associated data to file for restart. More...
 
void read (std::ifstream &restart_file)
 Read nodal position and associated data from file for restart. More...
 
virtual void pin_all ()
 The pin_all() function must be overloaded by SolidNodes, so we put the virtual interface here to avoid virtual functions in Data. More...
 
virtual void unpin_all ()
 The unpin_all() function must be overloaded by SolidNode, so we put the virtual interface here to avoid virtual functions in Data. More...
 
unsigned hang_code ()
 Code that encapsulates the hanging status of the node (incl. the geometric hanging status) as $ \sum_{i=-1}{nval-1} Node::is_hanging(i) 2^{i+1} $. More...
 
HangInfo *const & hanging_pt () const
 Return pointer to hanging node data (this refers to the geometric hanging node status) (const version). More...
 
HangInfo *const & hanging_pt (const int &i) const
 Return pointer to hanging node data for value i (const version) More...
 
bool is_hanging () const
 Test whether the node is geometrically hanging. More...
 
bool is_hanging (const int &i) const
 Test whether the i-th value is hanging. More...
 
void set_hanging_pt (HangInfo *const &hang_pt, const int &i)
 Set the hanging data for the i-th value. (hang_pt=0 to make non-hanging) More...
 
void set_nonhanging ()
 Label node as non-hanging node by removing all hanging node data. More...
 
void resize (const unsigned &n_value)
 Resize the number of equations. More...
 
virtual void constrain_positions ()
 Constrain the positions when the node is made hanging Empty virtual function that is overloaded in SolidNodes. More...
 
virtual void unconstrain_positions ()
 Unconstrain the positions when the node is made non-hanging Empty virtual function that is overloaded in SolidNodes. More...
 
virtual void make_periodic (Node *const &node_pt)
 Make the node periodic by copying the values from node_pt. Note that the coordinates will always remain independent, even though this may lead to (a little) unrequired information being stored. Broken virtual (only implemented in BoundaryNodes) More...
 
virtual void make_periodic_nodes (const Vector< Node * > &periodic_nodes_pt)
 Make the nodes passed in the vector periodic_nodes share the same data as this node. More...
 
virtual void get_boundaries_pt (std::set< unsigned > *&boundaries_pt)
 Return a pointer to set of mesh boundaries that this node occupies; this will be overloaded by BoundaryNodes. The default behaviour is that the Node does not lie on any boundaries so the pointer to the set of boundaries is NULL. More...
 
virtual bool is_on_boundary () const
 Test whether the Node lies on a boundary. The "bulk" Node cannot lie on a boundary, so return false. This will be overloaded by BoundaryNodes. More...
 
virtual bool is_on_boundary (const unsigned &b) const
 Test whether the node lies on mesh boundary b. The "bulk" Node cannot lie on a boundary, so return false. This will be overloaded by BoundaryNodes. More...
 
virtual void add_to_boundary (const unsigned &b)
 Broken interface for adding the node to the mesh boundary b Essentially here for error reporting. More...
 
virtual void remove_from_boundary (const unsigned &b)
 Broken interface for removing the node from the mesh boundary b Here to provide error reporting. More...
 
virtual unsigned ncoordinates_on_boundary (const unsigned &b)
 Get the number of boundary coordinates on mesh boundary b. Broken virtual interface provides run-time error checking. More...
 
virtual bool boundary_coordinates_have_been_set_up ()
 Have boundary coordinates been set up? Broken virtual interface provides run-time error checking. More...
 
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. Broken virtual interface provides run-time error checking. More...
 
virtual void set_coordinates_on_boundary (const unsigned &b, const unsigned &k, const Vector< double > &boundary_zeta)
 Set the vector of the k-th generalised boundary coordinates on mesh boundary b. Broken virtual interface provides run-time error checking. More...
 
virtual void get_coordinates_on_boundary (const unsigned &b, Vector< double > &boundary_zeta)
 Return the vector of coordinates on mesh boundary b Broken virtual interface provides run-time error checking. More...
 
virtual void set_coordinates_on_boundary (const unsigned &b, const Vector< double > &boundary_zeta)
 Set the vector of coordinates on mesh boundary b Broken virtual interface provides run-time error checking. More...
 
void set_obsolete ()
 Mark node as obsolete. More...
 
void set_non_obsolete ()
 Mark node as non-obsolete. More...
 
bool is_obsolete ()
 Test whether node is obsolete. More...
 
double raw_value (const unsigned &i) const
 Return the i-th value stored at the Node. This interface does NOT take the hanging status of the Node into account. More...
 
double raw_value (const unsigned &t, const unsigned &i) const
 Return the i-th value at time level t (t=0: present, t>0: previous). This interface does NOT take the hanging status of the Node into account. More...
 
double value (const unsigned &i) const
 Return i-th value (dofs or pinned) at this node either directly or via hanging node representation. Note that this REDFINES the interface in Data Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called. More...
 
double value (const unsigned &t, const unsigned &i) const
 Return i-th value at time level t (t=0: present, t>0: previous) either directly or via hanging node representation. Note that this REDEFINES the interface in Data Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called. More...
 
void value (Vector< double > &values) const
 Compute Vector of values for the Data value taking the hanging node status into account. Note that this REDEFINES the interface in Data Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called. More...
 
Vector< double > value () const
 Return vector of values calculated using value(vector). More...
 
void value (const unsigned &t, Vector< double > &values) const
 Compute Vector of values (dofs or pinned) in this data at time level t (t=0: present; t>0: previous). This interface explicitly takes the hanging status into account. Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called. More...
 
void position (Vector< double > &pos) const
 Compute Vector of nodal positions either directly or via hanging node representation. More...
 
Vector< double > position () const
 Return vector of position of node at current time. More...
 
void position (const unsigned &t, Vector< double > &pos) const
 Compute Vector of nodal position at time level t (t=0: current; t>0: previous timestep), either directly or via hanging node representation. More...
 
double position (const unsigned &i) const
 Return i-th nodal coordinate either directly or via hanging node representation. More...
 
double position (const unsigned &t, const unsigned &i) const
 Return i-th nodal coordinate at time level t (t=0: current; t>0: previous time level), either directly or via hanging node representation. More...
 
double position_gen (const unsigned &k, const unsigned &i) const
 Return generalised nodal coordinate either directly or via hanging node representation. More...
 
double position_gen (const unsigned &t, const unsigned &k, const unsigned &i) const
 Return generalised nodal coordinate at time level t (t=0: current; t>0: previous time level), either directly or via hanging node representation. More...
 
double dposition_dt (const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt, either directly or via hanging node representation. More...
 
double dposition_dt (const unsigned &j, const unsigned &i) const
 Return the i-th component of j-th derivative of nodal position: d^jx/dt^j either directly or via hanging node representation. More...
 
double dposition_gen_dt (const unsigned &k, const unsigned &i) const
 i-th component of time derivative (velocity) of the generalised position, dx(k,i)/dt. ‘Type’: k; Coordinate direction: i. This function uses the hanging node representation if necessary. More...
 
double dposition_gen_dt (const unsigned &j, const unsigned &k, const unsigned &i) const
 i-th component of j-th time derivative (velocity) of the generalised position, d^jx(k,i)/dt^j. ‘Type’: k; Coordinate direction: i. This function uses the hanging node representation if necessary More...
 
virtual void node_update (const bool &update_all_time_levels_for_new_node=false)
 Interface for functions that update the nodal position using algebraic remeshing strategies. The interface is common to SpineNodes, AlgebraicNodes and MacroElementNodeUpdateNodes. The default is that the node does not "update itself" i.e. it is fixed in space. When implemented, this function should also execute the Node's auxiliary node update function (if any). More...
 
void set_auxiliary_node_update_fct_pt (AuxNodeUpdateFctPt aux_node_update_fct_pt)
 Set pointer to auxiliary update function – this can be used to update any nodal values following the update of the nodal position. This is needed e.g. to update the no-slip condition on moving boundaries. More...
 
bool has_auxiliary_node_update_fct_pt ()
 Boolean to indicate if node has a pointer to and auxiliary update function. More...
 
void perform_auxiliary_node_update_fct ()
 Execute auxiliary update function (if any) – this can be used to update any nodal values following the update of the nodal position. This is needed e.g. to update the no-slip condition on moving boundaries. More...
 
virtual unsigned ngeom_data () const
 Return the number of geometric data that affect the nodal position. The default value is zero (node is stationary) More...
 
virtual Data ** all_geom_data_pt ()
 Return a pointer to an array of all (geometric) data that affect the nodal position. The default value is zero (node is stationary) More...
 
virtual unsigned ngeom_object () const
 Return the number of geometric objects that affect the nodal position. The default value is zero (node is stationary) More...
 
virtual GeomObject ** all_geom_object_pt ()
 Return a pointer to an array of all (geometric) objects that affect the nodal position. The default value is zero (node is stationary) More...
 
void output (std::ostream &outfile)
 Output nodal position. More...
 
void add_values_to_vector (Vector< double > &vector_of_values)
 Add all data and time history values to the vector. Overloaded to add the position information as well. More...
 
void read_values_from_vector (const Vector< double > &vector_of_values, unsigned &index)
 Read all data and time history values from the vector starting from index. On return the index will be set the value at the end of the data that has been read in Overload to also read the position information. More...
 
- Public Member Functions inherited from oomph::Data
virtual void clear_copied_pointers ()
 Helper function that should be overloaded derived classes that contain copies of data. The function must unset (NULL out) the internal pointers to the copied data. This is used when destructing data to ensure that all pointers remain valid. The default implementation throws an error because Data cannot be a copy. More...
 
 Data ()
 Default: Just set pointer to (steady) timestepper. No storage for values is allocated. More...
 
 Data (const unsigned &initial_n_value)
 Default constructor for steady problems: assign memory for initial_n_value values. More...
 
 Data (TimeStepper *const &time_stepper_pt, const unsigned &initial_n_value, const bool &allocate_storage=true)
 Constructor for unsteady problems: assign memory for initial_n_value values and any memory required by the Timestepper for the storage of history values. More...
 
 Data (const Data &data)=delete
 Broken copy constructor. More...
 
void operator= (const Data &)=delete
 Broken assignment operator. More...
 
virtual ~Data ()
 Destructor, deallocates memory assigned for data. More...
 
void set_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Set a new timestepper by resizing the appropriate storage. If already assigned the equation numbering will not be altered. More...
 
TimeStepper *& time_stepper_pt ()
 Return the pointer to the timestepper. More...
 
TimeStepper *const & time_stepper_pt () const
 Return the pointer to the timestepper (const version). More...
 
virtual bool is_a_copy () const
 Return a boolean to indicate whether the Data objact contains any copied values. A base Data object can never be a copy so the default implementation always returns false. More...
 
virtual bool is_a_copy (const unsigned &i) const
 Return flag to indicate whether the i-th value is a copy. A base Data object can never be a copy so the default implementation always returns false. More...
 
void set_value (const unsigned &i, const double &value_)
 Set the i-th stored data value to specified value. The only reason that we require an explicit set function is because we redefine value() in the Node class to interpolate the values for nodes that are hanging and so we cannot return a reference to the value in this case. More...
 
void set_value (const unsigned &t, const unsigned &i, const double &value_)
 Set the t-th history value of the i-th stored data value to specified value. More...
 
double value (const unsigned &i) const
 Return i-th stored value. This function is not virtual so that it can be inlined. This means that if we have an explicit pointer to a Data object Data* data_pt->value() always returns the "raw" stored value. More...
 
double value (const unsigned &t, const unsigned &i) const
 Return i-th value at time level t (t=0: present, t>0: previous) This function is not virtual so that it can be inlined. This means that if we have an explicit pointer to a Data object Data* data_pt->value() always returns to the "raw" stored value. More...
 
void value (Vector< double > &values) const
 Compute Vector of values for the Data value. More...
 
void value (const unsigned &t, Vector< double > &values) const
 Compute Vector of values (dofs or pinned) in this data at time level t (t=0: present; t>0: previous). More...
 
double * value_pt (const unsigned &i) const
 Return the pointer to the i-the stored value. Typically this is required when direct access to the stored value is required, e.g. when writing functions that return a reference to a variable that is stored in a Data object. More...
 
double * value_pt (const unsigned &t, const unsigned &i) const
 Return the pointer to the i-th stored value, or any of its history values (const version). Typically this is required when direct access to the stored value is required, e.g. when writing functions that return a reference to a variable that is stored in a Data object. More...
 
bool does_pointer_correspond_to_value (double *const &parameter_pt)
 Check whether the pointer parameter_pt addresses internal data values. More...
 
void copy (Data *orig_data_pt)
 Copy Data values from specified Data object. More...
 
void dump (std::ostream &dump_file) const
 Dump the data object to a file. More...
 
void read (std::ifstream &restart_file)
 Read data object from a file. More...
 
long * eqn_number_pt (const unsigned &i)
 Return the pointer to the equation number of the i-th stored variable. More...
 
long & eqn_number (const unsigned &i)
 Return the equation number of the i-th stored variable. More...
 
long eqn_number (const unsigned &i) const
 Return the equation number of the i-th stored variable. More...
 
void pin (const unsigned &i)
 Pin the i-th stored variable. More...
 
void unpin (const unsigned &i)
 Unpin the i-th stored variable. More...
 
void pin_all ()
 Pin all the stored variables. More...
 
void unpin_all ()
 Unpin all the stored variables. More...
 
bool is_pinned (const unsigned &i) const
 Test whether the i-th variable is pinned (1: true; 0: false). More...
 
bool is_segregated_solve_pinned (const unsigned &i)
 Test whether the i-th variable is temporaily pinned for a segregated solve. More...
 
void constrain (const unsigned &i)
 Constrain the i-th stored variable when making hanging data If the data is already pinned leave it along, otherwise mark as constrained (hanging) More...
 
void unconstrain (const unsigned &i)
 Unconstrain the i-th stored variable when make the data nonhanging. Only unconstrain if it was actually constrained (hanging) More...
 
void constrain_all ()
 Constrain all the stored variables when the data is made hanging. More...
 
void unconstrain_all ()
 Unconstrain all the stored variables when the data is made nonhanging. More...
 
bool is_constrained (const unsigned &i)
 Test whether the i-th variable is constrained (1: true; 0: false). More...
 
unsigned self_test ()
 Self-test: Have all values been classified as pinned/unpinned? Return 0 if OK. More...
 
unsigned nvalue () const
 Return number of values stored in data object (incl pinned ones). More...
 
unsigned ntstorage () const
 Return total number of doubles stored per value to record time history of each value (one for steady problems). More...
 
virtual void describe_dofs (std::ostream &out, const std::string &current_string) const
 Function to describe the dofs of the Node. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...) More...
 
virtual void add_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 Add pointers to all unpinned and unconstrained data to a map indexed by (global) equation number. More...
 
void set_halo (const unsigned &non_halo_proc_ID)
 Label the node as halo and specify processor that holds non-halo counterpart. More...
 
void set_nonhalo ()
 Label the node as not being a halo. More...
 
bool is_halo () const
 Is this Data a halo? More...
 
int non_halo_proc_ID ()
 ID of processor ID that holds non-halo counterpart of halo node; negative if not a halo. More...
 
virtual void add_eqn_numbers_to_vector (Vector< long > &vector_of_eqn_numbers)
 Add all equation numbers to the vector in the internal storage order. More...
 
virtual void read_eqn_numbers_from_vector (const Vector< long > &vector_of_eqn_numbers, unsigned &index)
 Read all equation numbers from the vector starting from index. On return the index will be set to the value at the end of the data that has been read in. More...
 

Static Public Attributes

static unsigned No_independent_position = 10
 Static "Magic number" used to indicate that there is no independent position in a periodic node. More...
 
- Static Public Attributes inherited from oomph::Data
static long Is_pinned = -1
 Static "Magic number" used in place of the equation number to indicate that the value is pinned. More...
 
static long Is_segregated_solve_pinned = -3
 Static "Magic number" used in place of the equation number to indicate that the value is pinned, but only for the duration of a segregated solve. More...
 
static long Is_unclassified = -10
 Static "Magic number" used in place of the equation number to denote a value that hasn't been classified as pinned or free. More...
 
static long Is_constrained = -2
 Static "Magic number" used in place of the equation number to indicate that the value is constrained because it is associated with non-conforming element boundaries — a hanging node — (and is therefore pinned) More...
 

Protected Member Functions

void x_gen_range_check (const unsigned &t, const unsigned &k, const unsigned &i) const
 Private function to check that the arguemnts to the position functions are in range. More...
 
double * x_position_pt (const unsigned &i)
 Direct access to the pointer to the i-th stored coordinate data. More...
 
- Protected Member Functions inherited from oomph::Data
virtual void reset_copied_pointers ()
 Helper function that should be overloaded in derived classes that can contain copies of Data. The function must reset the internal pointers to the copied data. This is used when resizing data to ensure that all the pointers remain valid. The default implementation throws an error beacause Data cannot be a copy. More...
 

Protected Attributes

double ** X_position
 Array of pointers to the data holding the Eulerian positions. The storage format must be the same as the internal data storage so that we can implement the functions x() in generality here without the need for virtual functions. The first index will be a flat array of position types and coordinates and the second will be the number of time history values at each position type. More...
 
TimeStepperPosition_time_stepper_pt
 Pointer to the timestepper associated with the position data. More...
 
HangInfo ** Hanging_pt
 C-style array of pointers to hanging node info. It's set to NULL if the node isn't hanging. The first entry (0) is the geometric hanging node data. The remaining entries correspond to the hanging data for the other values stored at the node. Usually, these entries will be the same as the geometric hanging node data represented by Hanging_pt[0], but this is not necessarily the case; e.g. the pressure in Taylor Hood has different hanging node data from the velocities. More...
 
unsigned Ndim
 Eulerian dimension of the node. More...
 
unsigned Nposition_type
 Number of coordinate types used in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements, etc). More...
 
bool Obsolete
 Flag to indicate that the Node has become obsolete — usually during mesh refinement process. More...
 
AuxNodeUpdateFctPt Aux_node_update_fct_pt
 Pointer to auxiliary update function – this can be used to update any nodal values following the update of the nodal position. This is needed e.g. to update the no-slip condition on moving boundaries. More...
 
- Protected Attributes inherited from oomph::Data
Data ** Copy_of_data_pt
 C-style array of any Data objects that contain copies of the current Data object's data values. More...
 
unsigned Ncopies
 Number of Data that contain copies of this Data object's values. More...
 

Friends

class BoundaryNodeBase
 to construct periodic Nodes More...
 
std::ostream & operator<< (std::ostream &out, const Node &d)
 Output operator: output location and all values at all times, along with any extra information stored for the timestepper. More...
 

Additional Inherited Members

- Static Protected Attributes inherited from oomph::Data
static TimeStepperDefault_static_time_stepper_pt = new Steady<0>()
 Default (static) timestepper used in steady problems. More...
 

Detailed Description

Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a given dimension.

The nodal coordinates are used in the elements' mapping between local and global coordinates and in the simplest case (stationary nodes in Lagrange-type elements) this mapping is given by

\[ x_i = \sum_{j=1}^{N_{node}} X_{ij} \psi_{j}(s_k) \]

so we need only access to the nodal coordinates $ X_{ij}\ (i=1..DIM) $ of all nodes $ j $ : provided by the Node member function

cstr elem_len * i
Definition: cfortran.h:603
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060

If the nodal positions are time-dependent, the mapping becomes

\[ x_i(t) = \sum_{j=1}^{N_{node}} X_{ij}(t) \ \psi_{j}(s_k). \]

Within the computation (where time is only evaluated at discrete levels) this becomes

\[ x_{ti} = \sum_{j=1}^{N_{node}} X_{ijt} \ \psi_{j}(s_k). \]

and we need access to the nodal coordinates $ X_{ijt} \ (i=1..DIM) $ of all nodes $ j $ at the present (t=0) and previous (t>0) timesteps: provided by the Node member function

char t
Definition: cfortran.h:568

Note: The interpretation of the history values is slightly more subtle than that. Depending on the positional TimeStepper used, only a limited number of the positional history values accessed Node::x(t,i) represent previous nodal positions; the others are generalised history values that the TimeStepper uses to determine approximations for the time-derivatives of the nodal positions.

Finally, some elements employ mappings that involve additional, generalised coordinates. For instance, in Hermite elements the mapping between local and global coordinates is based on an independent interpolation for the global coordinates and their derivative w.r.t. to the local coordinates. In such elements, the mapping becomes

\[ x_i = \sum_{j=1}^{N_{node}} \sum_{k=1}^{N_{type}} X_{ijk} \psi_{jk}(s_k) \]

where $ N_{type} $ is the number of the different types of generalised coordinates involved in the mapping. For instance, in 1D Hermite elements $ N_{type}=2 $ and k=0 corresponds to the global coordinates while k=1 corresponds to the derivative of the global coordinates w.r.t. to the local coordinate. In such cases we need access to the generalised nodal coordinates $ X_{ijk} \ (i=1..DIM, \ k=1..N_{type}) $ of all nodes $ j $. Access is provided by the Node member function

double & x_gen(const unsigned &k, const unsigned &i)
Reference to the generalised position x(k,i). ‘Type’: k; Coordinate direction: i.
Definition: nodes.h:1126

and the corresponding time-dependent version

While this is all pretty straightforward, it does make the argument list of the Node constructors rather lengthy.

Definition at line 905 of file nodes.h.

Member Typedef Documentation

◆ AuxNodeUpdateFctPt

typedef void(* oomph::Node::AuxNodeUpdateFctPt) (Node *)

Function pointer to auxiliary node update function.

Definition at line 909 of file nodes.h.

Constructor & Destructor Documentation

◆ Node() [1/4]

oomph::Node::Node ( )

Default constructor.

Definition at line 1574 of file nodes.cc.

References oomph::LeakCheckNames::Node_build.

◆ Node() [2/4]

oomph::Node::Node ( const unsigned &  n_dim,
const unsigned &  n_position_type,
const unsigned &  initial_n_value,
const bool &  allocate_x_position = true 
)

Steady constructor, for a Node of spatial dimension n_dim. Allocates storage for initial_n_value values. NPosition_type is the number of coordinate types needed in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements, etc).

Steady Constructor, allocates storage for initial_n_value values at a node of spatial dimension NDim. nposition_type: # of coordinate types needed in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements, etc).

Definition at line 1595 of file nodes.cc.

References oomph::LeakCheckNames::Node_build, and X_position.

◆ Node() [3/4]

oomph::Node::Node ( TimeStepper *const &  time_stepper_pt,
const unsigned &  n_dim,
const unsigned &  n_position_type,
const unsigned &  initial_n_value,
const bool &  allocate_x_position = true 
)

Unsteady constructor for a node of spatial dimension n_dim. Allocates storage for initial_n_value values with history values as required by the timestepper. n_position_type: # of coordinate types needed in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements).

Unsteady Constructor for a node of spatial dimension n_dim. Allocates storage for initial_n_value values with history values as required by timestepper. n_position_type: # of coordinate types needed in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements)

Definition at line 1644 of file nodes.cc.

References oomph::LeakCheckNames::Node_build, oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, t, and X_position.

◆ ~Node()

oomph::Node::~Node ( )
virtual

Destructor: Clean up the memory allocated for nodal position.

Destructor to clean up the memory allocated for nodal position.

Definition at line 1695 of file nodes.cc.

References Hanging_pt, oomph::LeakCheckNames::Node_build, oomph::Data::nvalue(), and X_position.

◆ Node() [4/4]

oomph::Node::Node ( const Node node)
delete

Broken copy constructor.

Member Function Documentation

◆ add_to_boundary()

void oomph::Node::add_to_boundary ( const unsigned &  b)
virtual

Broken interface for adding the node to the mesh boundary b Essentially here for error reporting.

Interface for function to add the node to the mesh boundary b. Broken here in order to report run-time errors. Must be overloaded by all boundary nodes.

Definition at line 2336 of file nodes.cc.

Referenced by oomph::Mesh::add_boundary_node(), and oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection().

◆ add_values_to_vector()

void oomph::Node::add_values_to_vector ( Vector< double > &  vector_of_values)
virtual

Add all data and time history values to the vector. Overloaded to add the position information as well.

Add position data and time-history values to the vector after the addition of the "standard" data stored at the node.

Reimplemented from oomph::Data.

Reimplemented in oomph::SolidNode.

Definition at line 2777 of file nodes.cc.

References oomph::Data::add_values_to_vector(), i, ndim(), nposition_type(), oomph::TimeStepper::ntstorage(), position_time_stepper_pt(), t, and X_position.

Referenced by oomph::SolidNode::add_values_to_vector(), oomph::Problem::get_data_to_be_sent_during_load_balancing(), and oomph::Problem::synchronise_dofs().

◆ all_geom_data_pt()

virtual Data** oomph::Node::all_geom_data_pt ( )
inlinevirtual

Return a pointer to an array of all (geometric) data that affect the nodal position. The default value is zero (node is stationary)

Reimplemented in oomph::SpineNode.

Definition at line 1632 of file nodes.h.

Referenced by oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data().

◆ all_geom_object_pt()

virtual GeomObject** oomph::Node::all_geom_object_pt ( )
inlinevirtual

Return a pointer to an array of all (geometric) objects that affect the nodal position. The default value is zero (node is stationary)

Reimplemented in oomph::SpineNode, oomph::MacroElementNodeUpdateNode, and oomph::AlgebraicNode.

Definition at line 1647 of file nodes.h.

Referenced by oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data().

◆ assign_eqn_numbers()

void oomph::Node::assign_eqn_numbers ( unsigned long &  global_ndof,
Vector< double * > &  dof_pt 
)
virtual

Assign global equation numbers; increment global number of unknowns, global_ndof; and add any new dofs to the dof_pt.

Assign (global) equation number. Overloaded version for nodes. Checks if a hanging value has a non-negative equation number and if so shouts and then sets it to Is_constrained. Then drops down to Data version which does the actual work.

Reimplemented from oomph::Data.

Reimplemented in oomph::SolidNode.

Definition at line 534 of file nodes.cc.

References oomph::Data::assign_eqn_numbers(), oomph::Data::constrain(), i, oomph::Data::is_constrained(), is_hanging(), oomph::Data::is_pinned(), and oomph::Data::nvalue().

◆ boundary_coordinates_have_been_set_up()

bool oomph::Node::boundary_coordinates_have_been_set_up ( )
virtual

Have boundary coordinates been set up? Broken virtual interface provides run-time error checking.

Interface for function to report if boundary coordinates have been set up for this node Broken here in order to report run-time errors. Must be overloaded by all boundary nodes.

Definition at line 2323 of file nodes.cc.

◆ constrain_positions()

virtual void oomph::Node::constrain_positions ( )
inlinevirtual

Constrain the positions when the node is made hanging Empty virtual function that is overloaded in SolidNodes.

Reimplemented in oomph::SolidNode.

Definition at line 1345 of file nodes.h.

Referenced by set_hanging_pt().

◆ copied_node_pt()

virtual Node* oomph::Node::copied_node_pt ( ) const
inlinevirtual

Return pointer to copied node (null if the current node is not a copy – always the case here; it's overloaded for boundary nodes)

Definition at line 1107 of file nodes.h.

Referenced by oomph::FourierDecomposedHelmholtzDtNBoundaryElement< ELEMENT >::complete_setup_of_dependencies(), oomph::HelmholtzDtNBoundaryElement< ELEMENT >::complete_setup_of_dependencies(), and oomph::BoundaryNodeBase::make_node_periodic().

◆ copy()

void oomph::Node::copy ( Node orig_node_pt)

◆ does_pointer_correspond_to_position_data()

virtual bool oomph::Node::does_pointer_correspond_to_position_data ( double *const &  parameter_pt)
inlinevirtual

Check whether the pointer parameter_pt addresses position data values. It never does for a standard node, because the positions are not data.

Reimplemented in oomph::SolidNode.

Definition at line 1042 of file nodes.h.

◆ dposition_dt() [1/2]

double oomph::Node::dposition_dt ( const unsigned &  i) const

Return the i-th component of nodal velocity: dx/dt, either directly or via hanging node representation.

Return the i-th component of nodal velocity: dx/dt / Use the hanging node representation if required.

Definition at line 2659 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), position(), Position_time_stepper_pt, t, and oomph::TimeStepper::weight().

Referenced by oomph::FSI_functions::apply_no_slip_on_moving_wall(), and oomph::FiniteElement::dnodal_position_dt().

◆ dposition_dt() [2/2]

double oomph::Node::dposition_dt ( const unsigned &  j,
const unsigned &  i 
) const

Return the i-th component of j-th derivative of nodal position: d^jx/dt^j either directly or via hanging node representation.

Return the i-th component of j-th derivative of nodal position: d^jx/dt^j. Use the hanging node representation.

Definition at line 2683 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), position(), Position_time_stepper_pt, t, and oomph::TimeStepper::weight().

◆ dposition_gen_dt() [1/2]

double oomph::Node::dposition_gen_dt ( const unsigned &  j,
const unsigned &  k,
const unsigned &  i 
) const

i-th component of j-th time derivative (velocity) of the generalised position, d^jx(k,i)/dt^j. ‘Type’: k; Coordinate direction: i. This function uses the hanging node representation if necessary

i-th component of j-th time derivative (velocity) of the generalised position, d^jx(k,i)/dt^j. ‘Type’: k; Coordinate direction: i. Use the hanging node representation.

Definition at line 2734 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), position_gen(), Position_time_stepper_pt, t, and oomph::TimeStepper::weight().

◆ dposition_gen_dt() [2/2]

double oomph::Node::dposition_gen_dt ( const unsigned &  k,
const unsigned &  i 
) const

i-th component of time derivative (velocity) of the generalised position, dx(k,i)/dt. ‘Type’: k; Coordinate direction: i. This function uses the hanging node representation if necessary.

i-th component of time derivative (velocity) of the generalised position, dx(k,i)/dt. ‘Type’: k; Coordinate direction: i. Use the hanging node representation

Definition at line 2708 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), position_gen(), Position_time_stepper_pt, t, and oomph::TimeStepper::weight().

Referenced by oomph::FiniteElement::dnodal_position_gen_dt().

◆ dump()

void oomph::Node::dump ( std::ostream &  dump_file) const
virtual

Dump nodal position and associated data to file for restart.

Dump nodal positions and associated data to file for restart.

Reimplemented in oomph::SolidNode.

Definition at line 1969 of file nodes.cc.

References oomph::Data::dump(), Ndim, Nposition_type, oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, t, and X_position.

Referenced by oomph::SolidNode::dump().

◆ dx_dt() [1/2]

double oomph::Node::dx_dt ( const unsigned &  i) const

Return the i-th component of nodal velocity: dx/dt.

Definition at line 1817 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, t, oomph::TimeStepper::weight(), and x().

Referenced by oomph::FiniteElement::raw_dnodal_position_dt().

◆ dx_dt() [2/2]

double oomph::Node::dx_dt ( const unsigned &  j,
const unsigned &  i 
) const

Return the i-th component of j-th derivative of nodal position: d^jx/dt^j.

Definition at line 1841 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, t, oomph::TimeStepper::weight(), and x().

◆ dx_gen_dt() [1/2]

double oomph::Node::dx_gen_dt ( const unsigned &  j,
const unsigned &  k,
const unsigned &  i 
) const

i-th component of j-th time derivative (velocity) of the generalised position, d^jx(k,i)/dt^j. ‘Type’: k; Coordinate direction: i.

Definition at line 1890 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, t, oomph::TimeStepper::weight(), and x_gen().

◆ dx_gen_dt() [2/2]

double oomph::Node::dx_gen_dt ( const unsigned &  k,
const unsigned &  i 
) const

i-th component of time derivative (velocity) of the generalised position, dx(k,i)/dt. ‘Type’: k; Coordinate direction: i.

Definition at line 1865 of file nodes.cc.

References i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, t, oomph::TimeStepper::weight(), and x_gen().

Referenced by oomph::FiniteElement::raw_dnodal_position_gen_dt().

◆ get_boundaries_pt()

virtual void oomph::Node::get_boundaries_pt ( std::set< unsigned > *&  boundaries_pt)
inlinevirtual

◆ get_coordinates_on_boundary() [1/2]

void oomph::Node::get_coordinates_on_boundary ( const unsigned &  b,
const unsigned &  k,
Vector< double > &  boundary_zeta 
)
virtual

Return the vector of the k-th generalised boundary coordinates on mesh boundary b. Broken virtual interface provides run-time error checking.

Interface for function to get the k-th generalised boundary coordinate of the node on boundary b. Broken here in order to provide run-time error reporting. Must be overloaded by all boundary nodes.

Definition at line 2379 of file nodes.cc.

Referenced by oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::TriangleMesh< ELEMENT >::compute_boundary_segments_connectivity_and_initial_zeta_values(), oomph::TriangleMeshBase::dump_triangulateio(), get_coordinates_on_boundary(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_load_balance_helper(), oomph::TriangleMesh< ELEMENT >::identify_boundary_segments_and_assign_initial_zeta_values(), oomph::TriangleMesh< ELEMENT >::output_boundary_coordinates(), oomph::TriangleMesh< ELEMENT >::re_assign_initial_zeta_values_for_internal_boundary(), oomph::TriangleMesh< ELEMENT >::re_scale_re_assigned_initial_zeta_values_for_internal_boundary(), oomph::RefineableTriangleMesh< ELEMENT >::send_boundary_node_info_of_shared_nodes(), oomph::UnstructuredTwoDMeshGeometryBase::setup_boundary_coordinates(), oomph::RefineableTetgenMesh< ELEMENT >::snap_nodes_onto_boundary(), oomph::RefineableTriangleMesh< ELEMENT >::snap_nodes_onto_boundary(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), oomph::UnstructuredTwoDMeshGeometryBase::snap_nodes_onto_geometric_objects(), oomph::TetMeshBase::snap_to_quadratic_surface(), oomph::TriangleMesh< ELEMENT >::synchronize_boundary_coordinates(), oomph::RefineableTriangleMesh< ELEMENT >::update_open_curve_after_restart(), oomph::RefineableTriangleMesh< ELEMENT >::update_open_curve_using_elements_area(), oomph::RefineableTriangleMesh< ELEMENT >::update_open_curve_using_face_mesh(), oomph::RefineableTriangleMesh< ELEMENT >::update_polygon_after_restart(), oomph::RefineableTriangleMesh< ELEMENT >::update_polygon_using_elements_area(), oomph::RefineableTriangleMesh< ELEMENT >::update_polygon_using_face_mesh(), oomph::FaceElement::zeta_nodal(), and oomph::GenericLagrangeInterpolatedProjectableElement< ELEMENT >::zeta_nodal().

◆ get_coordinates_on_boundary() [2/2]

virtual void oomph::Node::get_coordinates_on_boundary ( const unsigned &  b,
Vector< double > &  boundary_zeta 
)
inlinevirtual

Return the vector of coordinates on mesh boundary b Broken virtual interface provides run-time error checking.

Definition at line 1420 of file nodes.h.

References get_coordinates_on_boundary().

◆ hang_code()

unsigned oomph::Node::hang_code ( )
inline

Code that encapsulates the hanging status of the node (incl. the geometric hanging status) as $ \sum_{i=-1}{nval-1} Node::is_hanging(i) 2^{i+1} $.

Definition at line 1213 of file nodes.h.

References i, is_hanging(), and oomph::Data::nvalue().

Referenced by oomph::Mesh::doc_mesh_distribution().

◆ hanging_pt() [1/2]

HangInfo* const& oomph::Node::hanging_pt ( ) const
inline

Return pointer to hanging node data (this refers to the geometric hanging node status) (const version).

Definition at line 1228 of file nodes.h.

References Hanging_pt.

Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::TreeBasedRefineableMesh< ELEMENT >::additional_synchronise_hanging_nodes(), oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), oomph::RefineableSolidElement::assign_solid_hanging_local_eqn_numbers(), oomph::Mesh::check_halo_schemes(), oomph::TreeBasedRefineableMeshBase::complete_hanging_nodes(), oomph::TreeBasedRefineableMeshBase::complete_hanging_nodes_recursively(), oomph::Mesh::delete_all_external_storage(), oomph::RefineableAdvectionDiffusionEquations< DIM >::dinterpolated_u_adv_diff_ddata(), oomph::RefineableAxisymAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::RefineableSphericalAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::RefineableAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::RefineableGeneralisedNewtonianNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableSpaceTimeNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableSpaceTimeNavierStokesMixedOrderEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableYoungLaplaceEquations::fill_in_contribution_to_residuals(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePolarNavierStokesEquations::fill_in_generic_residual_contribution(), oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_cons_axisym_adv_diff(), oomph::RefineableLinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_linearised_axi_nst(), oomph::RefineableLinearisedNavierStokesEquations::fill_in_generic_residual_contribution_linearised_nst(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), oomph::RefineableSphericalAdvectionDiffusionEquations::fill_in_generic_residual_contribution_spherical_adv_diff(), oomph::RefineableSphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineableNavierStokesBoussinesqElement< NST_ELEMENT, AD_ELEMENT >::fill_in_off_diagonal_block_analytic(), oomph::RefineableAdvectionDiffusionBoussinesqElement< AD_ELEMENT, NST_ELEMENT >::fill_in_off_diagonal_block_analytic(), oomph::RefineableBuoyantQCrouzeixRaviartElement< DIM >::fill_in_off_diagonal_jacobian_blocks_analytic(), oomph::RefineablePseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd(), oomph::RefineableSolidElement::geom_data_pt(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineablePVDEquationsWithPressure< DIM >::get_mass_matrix_diagonal(), oomph::RefineableElement::identify_field_data_for_interactions(), oomph::RefineableSolidElement::identify_geometric_data(), oomph::RefineableGeneralisedNewtonianQTaylorHoodElement< DIM >::identify_load_data(), oomph::RefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::identify_load_data(), oomph::RefineableQTaylorHoodElement< DIM >::identify_load_data(), oomph::RefineableQCrouzeixRaviartElement< DIM >::identify_load_data(), oomph::RefineableQTaylorHoodSpaceTimeElement< DIM >::identify_load_data(), oomph::RefineableQTaylorHoodMixedOrderSpaceTimeElement< DIM >::identify_load_data(), oomph::RefineablePolarTaylorHoodElement::insert_load_data(), oomph::RefineablePolarCrouzeixRaviartElement::insert_load_data(), oomph::SolidNode::lagrangian_position(), oomph::SolidNode::lagrangian_position_gen(), oomph::RefineableSolidElement::ngeom_data(), oomph::AlgebraicMesh::node_update(), oomph::Mesh::node_update(), oomph::AlgebraicNode::node_update(), oomph::MacroElementNodeUpdateNode::node_update(), oomph::RefineableQElement< 3 >::oc_hang_helper(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::oc_hang_helper(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::MGSolver< DIM >::plot(), position(), position_gen(), oomph::RefineableQElement< 2 >::quad_hang_helper(), oomph::Missing_masters_functions::recursively_add_masters_of_external_haloed_node(), oomph::Multi_domain_functions::recursively_add_masters_of_external_haloed_node(), oomph::RefineableSolidTractionElement< ELEMENT >::refineable_fill_in_contribution_to_residuals_solid_traction(), oomph::RefineableImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::RefineableFSIImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(), oomph::RefineableNavierStokesFluxControlElement< ELEMENT >::refineable_fill_in_generic_residual_contribution_fluid_traction(), oomph::Problem::remove_duplicate_data(), oomph::MGSolver< DIM >::setup_interpolation_matrices_unstructured(), oomph::HelmholtzMGPreconditioner< DIM >::setup_interpolation_matrices_unstructured(), oomph::TreeBasedRefineableMeshBase::synchronise_hanging_nodes(), and value().

◆ hanging_pt() [2/2]

HangInfo* const& oomph::Node::hanging_pt ( const int &  i) const
inline

Return pointer to hanging node data for value i (const version)

Definition at line 1243 of file nodes.h.

References oomph::MPI_Helpers::communicator_pt(), Hanging_pt, i, ndim(), oomph::Data::nvalue(), and x().

◆ has_auxiliary_node_update_fct_pt()

bool oomph::Node::has_auxiliary_node_update_fct_pt ( )
inline

◆ is_hanging() [1/2]

bool oomph::Node::is_hanging ( ) const
inline

Test whether the node is geometrically hanging.

Definition at line 1285 of file nodes.h.

References Hanging_pt.

Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::TreeBasedRefineableMesh< ELEMENT >::additional_synchronise_hanging_nodes(), oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data(), assign_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), oomph::ExtrudedCubeMeshFromQuadMesh< ELEMENT >::build_mesh(), oomph::Mesh::check_halo_schemes(), oomph::TreeBasedRefineableMeshBase::complete_hanging_nodes(), oomph::TreeBasedRefineableMeshBase::complete_hanging_nodes_recursively(), oomph::Mesh::delete_all_external_storage(), oomph::RefineableAdvectionDiffusionEquations< DIM >::dinterpolated_u_adv_diff_ddata(), oomph::RefineableAxisymAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::RefineableSphericalAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::RefineableAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::RefineableGeneralisedNewtonianNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableSpaceTimeNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableSpaceTimeNavierStokesMixedOrderEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::Mesh::doc_shared_nodes(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePolarNavierStokesEquations::fill_in_generic_residual_contribution(), oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_cons_axisym_adv_diff(), oomph::RefineableLinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_linearised_axi_nst(), oomph::RefineableLinearisedNavierStokesEquations::fill_in_generic_residual_contribution_linearised_nst(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), oomph::RefineableSphericalAdvectionDiffusionEquations::fill_in_generic_residual_contribution_spherical_adv_diff(), oomph::RefineableSphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineableNavierStokesBoussinesqElement< NST_ELEMENT, AD_ELEMENT >::fill_in_off_diagonal_block_analytic(), oomph::RefineableAdvectionDiffusionBoussinesqElement< AD_ELEMENT, NST_ELEMENT >::fill_in_off_diagonal_block_analytic(), oomph::RefineableBuoyantQCrouzeixRaviartElement< DIM >::fill_in_off_diagonal_jacobian_blocks_analytic(), oomph::RefineablePseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineablePVDEquationsWithPressure< DIM >::get_mass_matrix_diagonal(), hang_code(), oomph::RefineableElement::identify_field_data_for_interactions(), oomph::RefineableGeneralisedNewtonianQTaylorHoodElement< DIM >::identify_load_data(), oomph::RefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::identify_load_data(), oomph::RefineableQTaylorHoodElement< DIM >::identify_load_data(), oomph::RefineableQCrouzeixRaviartElement< DIM >::identify_load_data(), oomph::RefineableQTaylorHoodSpaceTimeElement< DIM >::identify_load_data(), oomph::RefineableQTaylorHoodMixedOrderSpaceTimeElement< DIM >::identify_load_data(), oomph::RefineablePolarTaylorHoodElement::insert_load_data(), oomph::RefineablePolarCrouzeixRaviartElement::insert_load_data(), is_hanging(), oomph::SolidNode::lagrangian_position(), oomph::SolidNode::lagrangian_position_gen(), oomph::AlgebraicMesh::node_update(), oomph::Mesh::node_update(), oomph::AlgebraicNode::node_update(), oomph::MacroElementNodeUpdateNode::node_update(), oomph::RefineableQElement< 3 >::oc_hang_helper(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::RefineableAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableGeneralisedNewtonianAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableGeneralisedNewtonianQTaylorHoodElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableLinearisedAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableLinearisedQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineablePolarTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodSpaceTimeElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodMixedOrderSpaceTimeElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQSphericalTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::pin_elemental_redundant_nodal_solid_pressures(), oomph::MGSolver< DIM >::plot(), position(), position_gen(), oomph::RefineableQElement< 2 >::quad_hang_helper(), oomph::Missing_masters_functions::recursively_add_masters_of_external_haloed_node(), oomph::Multi_domain_functions::recursively_add_masters_of_external_haloed_node(), oomph::RefineableSolidTractionElement< ELEMENT >::refineable_fill_in_contribution_to_residuals_solid_traction(), oomph::RefineableImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::RefineableFSIImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(), oomph::RefineableNavierStokesFluxControlElement< ELEMENT >::refineable_fill_in_generic_residual_contribution_fluid_traction(), oomph::Problem::remove_duplicate_data(), oomph::MGSolver< DIM >::setup_interpolation_matrices_unstructured(), oomph::HelmholtzMGPreconditioner< DIM >::setup_interpolation_matrices_unstructured(), oomph::TreeBasedRefineableMeshBase::synchronise_hanging_nodes(), oomph::TreeBasedRefineableMeshBase::synchronise_nonhanging_nodes(), oomph::AxisymmetricTTaylorHoodElement::unpin_proper_nodal_pressure_dofs(), oomph::GeneralisedNewtonianAxisymmetricTTaylorHoodElement::unpin_proper_nodal_pressure_dofs(), oomph::GeneralisedNewtonianTTaylorHoodElement< DIM >::unpin_proper_nodal_pressure_dofs(), oomph::TTaylorHoodElement< DIM >::unpin_proper_nodal_pressure_dofs(), value(), and oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh().

◆ is_hanging() [2/2]

bool oomph::Node::is_hanging ( const int &  i) const
inline

Test whether the i-th value is hanging.

Definition at line 1298 of file nodes.h.

References Hanging_pt, i, is_hanging(), and oomph::Data::nvalue().

◆ is_obsolete()

bool oomph::Node::is_obsolete ( )
inline

Test whether node is obsolete.

Definition at line 1448 of file nodes.h.

References Obsolete.

Referenced by oomph::Mesh::prune_halo_elements_and_nodes().

◆ is_on_boundary() [1/2]

virtual bool oomph::Node::is_on_boundary ( ) const
inlinevirtual

Test whether the Node lies on a boundary. The "bulk" Node cannot lie on a boundary, so return false. This will be overloaded by BoundaryNodes.

Definition at line 1373 of file nodes.h.

Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::BackwardStepQuadMesh< ELEMENT >::build_mesh(), oomph::ExtrudedCubeMeshFromQuadMesh< ELEMENT >::build_mesh(), oomph::BrickFromTetMesh< ELEMENT >::build_mesh(), oomph::TriangleMesh< ELEMENT >::check_connections_of_polyline_nodes(), oomph::Mesh::copy_boundary_node_data_from_nodes(), oomph::Missing_masters_functions::get_required_master_nodal_information_helper(), oomph::Multi_domain_functions::get_required_master_nodal_information_helper(), oomph::Missing_masters_functions::get_required_nodal_information_helper(), oomph::Multi_domain_functions::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_load_balance_helper(), oomph::Edge::is_on_boundary(), oomph::TFace::is_on_boundary(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::MGSolver< DIM >::plot(), oomph::RefineableQuadMeshWithMovingCylinder< ELEMENT >::RefineableQuadMeshWithMovingCylinder(), oomph::RefineableTriangleMesh< ELEMENT >::send_boundary_node_info_of_shared_nodes(), oomph::RefineableTriangleMesh< ELEMENT >::snap_nodes_onto_boundary(), oomph::TetMeshBase::snap_to_quadratic_surface(), and oomph::TwoDAnnularMesh< ELEMENT >::wrap_into_annular_shape().

◆ is_on_boundary() [2/2]

virtual bool oomph::Node::is_on_boundary ( const unsigned &  b) const
inlinevirtual

Test whether the node lies on mesh boundary b. The "bulk" Node cannot lie on a boundary, so return false. This will be overloaded by BoundaryNodes.

Definition at line 1381 of file nodes.h.

◆ make_periodic()

void oomph::Node::make_periodic ( Node *const &  node_pt)
virtual

Make the node periodic by copying the values from node_pt. Note that the coordinates will always remain independent, even though this may lead to (a little) unrequired information being stored. Broken virtual (only implemented in BoundaryNodes)

Make the node periodic by copying values from node_pt. Broken virtual (only implemented in BoundaryNodes)

Definition at line 2257 of file nodes.cc.

Referenced by oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), and oomph::PeriodicOrbitTemporalMesh< ELEMENT >::PeriodicOrbitTemporalMesh().

◆ make_periodic_nodes()

void oomph::Node::make_periodic_nodes ( const Vector< Node * > &  periodic_nodes_pt)
virtual

Make the nodes passed in the vector periodic_nodes share the same data as this node.

Make the nodes passed in periodic_nodes_pt periodic by copying values across from this node. At present all the positions will be assumed to be independent. Broken virtual (only implemented in BoundaryNodes)

Definition at line 2270 of file nodes.cc.

◆ ncoordinates_on_boundary()

unsigned oomph::Node::ncoordinates_on_boundary ( const unsigned &  b)
virtual

Get the number of boundary coordinates on mesh boundary b. Broken virtual interface provides run-time error checking.

Interface to get the number of boundary coordinates on mesh boundary b. Broken here in order to provide run-time error reporting. Must be overloaded by all boundary nodes.

Definition at line 2363 of file nodes.cc.

Referenced by oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), and oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh().

◆ ndim()

unsigned oomph::Node::ndim ( ) const
inline

Return (Eulerian) spatial dimension of the node.

Definition at line 1054 of file nodes.h.

References Ndim.

Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::actions_before_newton_solve(), oomph::Problem::adapt(), oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::SurfactantTransportInterfaceElement::add_additional_residual_contributions_interface(), add_values_to_vector(), oomph::AdvectionDiffusionFluxElement< ELEMENT >::AdvectionDiffusionFluxElement(), oomph::FSI_functions::apply_no_slip_on_moving_wall(), oomph::Newmark< NSTEPS >::assign_initial_data_values(), oomph::Steady< NSTEPS >::assign_initial_positions_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_positions_impulsive(), oomph::BDF< NSTEPS >::assign_initial_positions_impulsive(), oomph::GeneralisedElement::assign_local_eqn_numbers(), oomph::Multi_domain_functions::aux_setup_multi_domain_interaction(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::backup(), oomph::SolidICProblem::backup_original_state(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::MeshAsGeomObject::build_it(), oomph::BDF< NSTEPS >::calculate_predicted_positions(), oomph::TR::calculate_predicted_positions(), oomph::Mesh::check_for_repeated_nodes(), oomph::Mesh::check_halo_schemes(), oomph::FiniteElement::check_jacobian(), oomph::Multi_domain_functions::construct_new_external_halo_master_node_helper(), oomph::Missing_masters_functions::construct_new_external_halo_master_node_helper(), oomph::Multi_domain_functions::construct_new_external_halo_node_helper(), oomph::Missing_masters_functions::construct_new_external_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_node_load_balance_helper(), copy(), oomph::Mesh::delete_all_external_storage(), oomph::Mesh::distribute(), oomph::Mesh::doc_boundary_coordinates(), oomph::FSI_functions::doc_fsi(), oomph::Mesh::doc_mesh_distribution(), oomph::Mesh::doc_shared_nodes(), oomph::DummyErrorEstimator::DummyErrorEstimator(), oomph::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data(), oomph::Problem::get_data_to_be_sent_during_load_balancing(), oomph::ElementWithMovingNodes::get_dnodal_coordinates_dgeom_dofs(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::RefineableElement::get_dresidual_dnodal_coordinates(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::get_field(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::get_field(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::get_field(), oomph::Missing_masters_functions::get_required_master_nodal_information_helper(), oomph::Multi_domain_functions::get_required_master_nodal_information_helper(), oomph::Missing_masters_functions::get_required_nodal_information_helper(), oomph::Multi_domain_functions::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_load_balance_helper(), oomph::FourierDecomposedHelmholtzBCElementBase< ELEMENT >::global_power_contribution(), oomph::PMLFourierDecomposedHelmholtzPowerMonitorElement< ELEMENT >::global_power_contribution(), hanging_pt(), oomph::HelmholtzBCElementBase< ELEMENT >::HelmholtzBCElementBase(), oomph::HelmholtzFluxElement< ELEMENT >::HelmholtzFluxElement(), oomph::HelmholtzFluxFromNormalDisplacementBCElement< HELMHOLTZ_BULK_ELEMENT, ELASTICITY_BULK_ELEMENT >::HelmholtzFluxFromNormalDisplacementBCElement(), oomph::LinearisedAxisymPoroelasticBJS_FSIElement< FLUID_BULK_ELEMENT, POROELASTICITY_BULK_ELEMENT >::LinearisedAxisymPoroelasticBJS_FSIElement(), oomph::LinearisedFSIAxisymmetricNStNoSlipBCElementElement< FLUID_BULK_ELEMENT, SOLID_BULK_ELEMENT >::LinearisedFSIAxisymmetricNStNoSlipBCElementElement(), oomph::LinearWaveFluxElement< ELEMENT >::LinearWaveFluxElement(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::local_equation(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::local_equation(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::local_equation(), oomph::NavierStokesFluxControlElement< ELEMENT >::NavierStokesFluxControlElement(), oomph::NavierStokesMixedOrderSpaceTimeTractionElement< ELEMENT >::NavierStokesMixedOrderSpaceTimeTractionElement(), oomph::NavierStokesSpaceTimeTractionElement< ELEMENT >::NavierStokesSpaceTimeTractionElement(), oomph::NavierStokesSurfaceDragTorqueElement< ELEMENT >::NavierStokesSurfaceDragTorqueElement(), oomph::NavierStokesSurfacePowerElement< ELEMENT >::NavierStokesSurfacePowerElement(), oomph::NavierStokesTractionElement< ELEMENT >::NavierStokesTractionElement(), oomph::NodeOrdering::node_global_position_comparison(), oomph::AlgebraicMesh::node_update(), oomph::Mesh::node_update(), oomph::MacroElementNodeUpdateNode::node_update(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()(), oomph::LinearElasticitySmoothMesh< LINEAR_ELASTICITY_ELEMENT >::operator()(), oomph::PoissonSmoothMesh< POISSON_ELEMENT >::operator()(), oomph::HermiteBeamElement::output(), oomph::DummyFaceElement< ELEMENT >::output(), output(), oomph::TriangleMesh< ELEMENT >::output_boundary_coordinates(), oomph::ElementWithExternalElement::output_external_elements(), oomph::FiniteElement::output_paraview(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), oomph::Tree::p_refine_if_required(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::pin_all(), oomph::MGSolver< DIM >::plot(), oomph::PMLHelmholtzFluxElement< ELEMENT >::PMLHelmholtzFluxElement(), oomph::PMLHelmholtzFluxFromNormalDisplacementBCElement< HELMHOLTZ_BULK_ELEMENT, ELASTICITY_BULK_ELEMENT >::PMLHelmholtzFluxFromNormalDisplacementBCElement(), oomph::PMLHelmholtzPowerElement< ELEMENT >::PMLHelmholtzPowerElement(), oomph::PoissonFluxElement< ELEMENT >::PoissonFluxElement(), oomph::PolarNavierStokesTractionElement< ELEMENT >::PolarNavierStokesTractionElement(), oomph::PolarStressIntegralElement< ELEMENT >::PolarStressIntegralElement(), position(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project(), oomph::Data::read(), read_values_from_vector(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::RefineableQSpectralElement< 3 >::rebuild_from_sons(), oomph::RefineableQSpectralElement< 2 >::rebuild_from_sons(), oomph::TreeBasedRefineableMeshBase::refine_as_in_reference_mesh(), oomph::Problem::remove_duplicate_data(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::reset(), oomph::SolidICProblem::reset_original_state(), oomph::Mesh::resize_halo_nodes(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::restore_positions(), oomph::Mesh::scale_mesh(), oomph::AlgebraicNode::self_test(), oomph::Problem::set_pinned_values_to_zero(), set_position_time_stepper(), oomph::MGSolver< DIM >::set_self_test_vector(), oomph::SolidICProblem::setup_problem(), oomph::IMRBase::shift_time_positions(), oomph::Steady< NSTEPS >::shift_time_positions(), oomph::Newmark< NSTEPS >::shift_time_positions(), oomph::NewmarkBDF< NSTEPS >::shift_time_positions(), oomph::BDF< NSTEPS >::shift_time_positions(), oomph::Tree::split_if_required(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::store_positions(), oomph::TetMeshVertex::TetMeshVertex(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::unpin_all(), oomph::UnsteadyHeatFluxElement< ELEMENT >::UnsteadyHeatFluxElement(), oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh(), and oomph::YoungLaplaceContactAngleElement< ELEMENT >::YoungLaplaceContactAngleElement().

◆ ngeom_data()

virtual unsigned oomph::Node::ngeom_data ( ) const
inlinevirtual

Return the number of geometric data that affect the nodal position. The default value is zero (node is stationary)

Reimplemented in oomph::SpineNode.

Definition at line 1625 of file nodes.h.

Referenced by oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data().

◆ ngeom_object()

virtual unsigned oomph::Node::ngeom_object ( ) const
inlinevirtual

Return the number of geometric objects that affect the nodal position. The default value is zero (node is stationary)

Reimplemented in oomph::SpineNode, oomph::MacroElementNodeUpdateNode, and oomph::AlgebraicNode.

Definition at line 1639 of file nodes.h.

Referenced by oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data().

◆ node_update()

virtual void oomph::Node::node_update ( const bool &  update_all_time_levels_for_new_node = false)
inlinevirtual

Interface for functions that update the nodal position using algebraic remeshing strategies. The interface is common to SpineNodes, AlgebraicNodes and MacroElementNodeUpdateNodes. The default is that the node does not "update itself" i.e. it is fixed in space. When implemented, this function should also execute the Node's auxiliary node update function (if any).

Reimplemented in oomph::SpineNode, oomph::SolidNode, oomph::MacroElementNodeUpdateNode, and oomph::AlgebraicNode.

Definition at line 1586 of file nodes.h.

Referenced by oomph::FiniteElement::node_update().

◆ nposition_type()

unsigned oomph::Node::nposition_type ( ) const
inline

Number of coordinate types needed in the mapping between local and global coordinates.

Definition at line 1016 of file nodes.h.

References Nposition_type.

Referenced by add_values_to_vector(), oomph::Steady< NSTEPS >::assign_initial_positions_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_positions_impulsive(), oomph::BDF< NSTEPS >::assign_initial_positions_impulsive(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), oomph::RefineableQElement< 2 >::build(), copy(), oomph::Missing_masters_functions::get_required_master_nodal_information_helper(), oomph::Multi_domain_functions::get_required_master_nodal_information_helper(), oomph::PRefineableQElement< 1, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::pin_all(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::pin_dofs_of_coordinate(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project(), read_values_from_vector(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::restore_positions(), oomph::Problem::set_pinned_values_to_zero(), set_position_time_stepper(), oomph::IMRBase::shift_time_positions(), oomph::Steady< NSTEPS >::shift_time_positions(), oomph::Newmark< NSTEPS >::shift_time_positions(), oomph::NewmarkBDF< NSTEPS >::shift_time_positions(), oomph::BDF< NSTEPS >::shift_time_positions(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::store_positions(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::unpin_all(), and oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::unpin_dofs_of_coordinate().

◆ operator=()

void oomph::Node::operator= ( const Node )
delete

Broken assignment operator.

◆ output()

void oomph::Node::output ( std::ostream &  outfile)

Output nodal position.

Output nodal coordinates.

Definition at line 2760 of file nodes.cc.

References i, ndim(), and x().

◆ perform_auxiliary_node_update_fct()

void oomph::Node::perform_auxiliary_node_update_fct ( )
inline

Execute auxiliary update function (if any) – this can be used to update any nodal values following the update of the nodal position. This is needed e.g. to update the no-slip condition on moving boundaries.

Definition at line 1615 of file nodes.h.

References Aux_node_update_fct_pt.

Referenced by oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineablePseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::RefineableElement::get_dresidual_dnodal_coordinates(), oomph::NavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::SpaceTimeNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableSpaceTimeNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::AlgebraicMesh::node_update(), and oomph::SolidNode::node_update().

◆ pin_all()

virtual void oomph::Node::pin_all ( )
inlinevirtual

The pin_all() function must be overloaded by SolidNodes, so we put the virtual interface here to avoid virtual functions in Data.

Reimplemented in oomph::SolidNode.

Definition at line 1197 of file nodes.h.

References oomph::Data::pin_all().

Referenced by oomph::SolidNode::pin_all().

◆ position() [1/5]

Vector<double> oomph::Node::position ( ) const
inline

Return vector of position of node at current time.

Definition at line 1525 of file nodes.h.

References ndim().

Referenced by dposition_dt(), and position().

◆ position() [2/5]

double oomph::Node::position ( const unsigned &  i) const

Return i-th nodal coordinate either directly or via hanging node representation.

Definition at line 2528 of file nodes.cc.

References hanging_pt(), i, is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), and x().

◆ position() [3/5]

double oomph::Node::position ( const unsigned &  t,
const unsigned &  i 
) const

Return i-th nodal coordinate at time level t (t=0: current; t>0: previous time level), either directly or via hanging node representation.

Definition at line 2560 of file nodes.cc.

References hanging_pt(), i, is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), t, and x().

◆ position() [4/5]

void oomph::Node::position ( const unsigned &  t,
Vector< double > &  pos 
) const

Compute Vector of nodal position at time level t (t=0: current; t>0: previous timestep), either directly or via hanging node representation.

Compute Vector of nodal position at timestep t (t=0: current; t>0: previous timestep), either directly or via hanging node representation.

Definition at line 2514 of file nodes.cc.

References i, ndim(), position(), and t.

◆ position() [5/5]

void oomph::Node::position ( Vector< double > &  pos) const

◆ position_gen() [1/2]

double oomph::Node::position_gen ( const unsigned &  k,
const unsigned &  i 
) const

Return generalised nodal coordinate either directly or via hanging node representation.

Definition at line 2592 of file nodes.cc.

References hanging_pt(), i, is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), and x_gen().

Referenced by dposition_gen_dt(), and oomph::FiniteElement::nodal_position_gen().

◆ position_gen() [2/2]

double oomph::Node::position_gen ( const unsigned &  t,
const unsigned &  k,
const unsigned &  i 
) const

Return generalised nodal coordinate at time level t (t=0: current; t>0: previous time level), either directly or via hanging node representation.

Definition at line 2624 of file nodes.cc.

References hanging_pt(), i, is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), t, and x_gen().

◆ position_is_a_copy() [1/2]

virtual bool oomph::Node::position_is_a_copy ( ) const
inlinevirtual

◆ position_is_a_copy() [2/2]

virtual bool oomph::Node::position_is_a_copy ( const unsigned &  i) const
inlinevirtual

Return whether the position coordinate i has been copied (always false)

Reimplemented in oomph::SolidNode.

Definition at line 1119 of file nodes.h.

◆ position_time_stepper_pt() [1/2]

TimeStepper*& oomph::Node::position_time_stepper_pt ( )
inline

Return a pointer to the position timestepper.

Definition at line 1022 of file nodes.h.

References Position_time_stepper_pt.

Referenced by add_values_to_vector(), oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::algebraic_node_update(), copy(), oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::ProjectableAdvectionDiffusionReactionElement< ADR_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymLinearElasticityElement< AXISYM_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricPoroelasticityElement< AXISYMMETRIC_POROELASTICITY_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableDarcyElement< DARCY_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableDisplacementBasedFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableFourierDecomposedHelmholtzElement< FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableGeneralisedNewtonianTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableGeneralisedNewtonianCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GenericLagrangeInterpolatedProjectableElement< ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableLinearElasticityElement< LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLFourierDecomposedHelmholtzElement< FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePoissonElement< POISSON_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePVDElement< PVD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePVDElementWithContinuousPressure< PVD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatMixedOrderSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_central_region(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_central_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_lower_right_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_upper_left_box(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_lower_right_region(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_upper_left_region(), oomph::ProjectableUnsteadyHeatSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::ProjectableUnsteadyHeatMixedOrderSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::output(), read_values_from_vector(), oomph::ContinuationStorageScheme::set_consistent_pinned_positions(), set_position_time_stepper(), oomph::SolidNode::set_position_time_stepper(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), and oomph::UnstructuredTwoDMeshGeometryBase::snap_nodes_onto_geometric_objects().

◆ position_time_stepper_pt() [2/2]

TimeStepper* const& oomph::Node::position_time_stepper_pt ( ) const
inline

Return a pointer to the position timestepper (const version).

Definition at line 1028 of file nodes.h.

References Position_time_stepper_pt.

◆ raw_value() [1/2]

double oomph::Node::raw_value ( const unsigned &  i) const
inline

Return the i-th value stored at the Node. This interface does NOT take the hanging status of the Node into account.

Definition at line 1455 of file nodes.h.

References i, and oomph::Data::value().

Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::FiniteElement::raw_nodal_value(), and value().

◆ raw_value() [2/2]

double oomph::Node::raw_value ( const unsigned &  t,
const unsigned &  i 
) const
inline

Return the i-th value at time level t (t=0: present, t>0: previous). This interface does NOT take the hanging status of the Node into account.

Definition at line 1463 of file nodes.h.

References i, t, and oomph::Data::value().

◆ read()

void oomph::Node::read ( std::ifstream &  restart_file)

Read nodal position and associated data from file for restart.

Read nodal positions and associated data from file for restart.

Definition at line 1996 of file nodes.cc.

References Ndim, Nposition_type, oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, oomph::Data::read(), oomph::Global_string_for_annotation::string(), t, and X_position.

Referenced by oomph::Mesh::read(), and oomph::SolidNode::read().

◆ read_values_from_vector()

void oomph::Node::read_values_from_vector ( const Vector< double > &  vector_of_values,
unsigned &  index 
)
virtual

Read all data and time history values from the vector starting from index. On return the index will be set the value at the end of the data that has been read in Overload to also read the position information.

Read the position data and its time histories from the vector after reading the "standard" data.

Reimplemented from oomph::Data.

Reimplemented in oomph::SolidNode.

Definition at line 2832 of file nodes.cc.

References i, ndim(), nposition_type(), oomph::TimeStepper::ntstorage(), position_time_stepper_pt(), oomph::Data::read_values_from_vector(), t, and X_position.

Referenced by oomph::SolidNode::read_values_from_vector(), oomph::Problem::send_data_to_be_sent_during_load_balancing(), and oomph::Problem::synchronise_dofs().

◆ remove_from_boundary()

void oomph::Node::remove_from_boundary ( const unsigned &  b)
virtual

Broken interface for removing the node from the mesh boundary b Here to provide error reporting.

Interface for function to remove the node from the mesh boundary b. Broken here in order to report run-time erorrs. Must be overloaded by all boundary nodes.

Definition at line 2350 of file nodes.cc.

Referenced by oomph::RefineableTriangleMesh< ELEMENT >::adapt(), oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::RefineableQuadMeshWithMovingCylinder< ELEMENT >::RefineableQuadMeshWithMovingCylinder(), oomph::Mesh::remove_boundary_node(), and oomph::Mesh::remove_boundary_nodes().

◆ resize()

void oomph::Node::resize ( const unsigned &  n_value)
virtual

◆ set_auxiliary_node_update_fct_pt()

void oomph::Node::set_auxiliary_node_update_fct_pt ( AuxNodeUpdateFctPt  aux_node_update_fct_pt)
inline

Set pointer to auxiliary update function – this can be used to update any nodal values following the update of the nodal position. This is needed e.g. to update the no-slip condition on moving boundaries.

Definition at line 1596 of file nodes.h.

References Aux_node_update_fct_pt.

◆ set_coordinates_on_boundary() [1/2]

void oomph::Node::set_coordinates_on_boundary ( const unsigned &  b,
const unsigned &  k,
const Vector< double > &  boundary_zeta 
)
virtual

Set the vector of the k-th generalised boundary coordinates on mesh boundary b. Broken virtual interface provides run-time error checking.

Interface for function to set the k-th generalised boundary coordinate of the node on boundary b. Broken here to provide run-time error reports. Must be overloaded by all boundary nodes.

Definition at line 2394 of file nodes.cc.

Referenced by oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::BrickFromTetMesh< ELEMENT >::build_mesh(), oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::ChannelWithLeafletMesh< ELEMENT >::ChannelWithLeafletMesh(), oomph::CollapsibleChannelMesh< ELEMENT >::CollapsibleChannelMesh(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_node_load_balance_helper(), oomph::CylinderWithFlagMesh< ELEMENT >::CylinderWithFlagMesh(), oomph::ElasticRefineableRectangularQuadMesh< ELEMENT >::ElasticRefineableRectangularQuadMesh(), oomph::FSIDrivenCavityMesh< ELEMENT >::FSIDrivenCavityMesh(), oomph::TriangleMesh< ELEMENT >::identify_boundary_segments_and_assign_initial_zeta_values(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), oomph::QuarterPipeMesh< ELEMENT >::QuarterPipeMesh(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::TriangleMeshBase::remesh_from_triangulateio(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::reposition_spines(), oomph::RefineableTriangleMesh< ELEMENT >::send_boundary_node_info_of_shared_nodes(), set_coordinates_on_boundary(), oomph::TetMeshBase::setup_boundary_coordinates(), oomph::XdaTetMesh< ELEMENT >::setup_boundary_coordinates(), oomph::UnstructuredTwoDMeshGeometryBase::setup_boundary_coordinates(), oomph::RefineableTriangleMesh< ELEMENT >::snap_nodes_onto_boundary(), oomph::TetMeshBase::snap_to_quadratic_surface(), oomph::TriangleMesh< ELEMENT >::synchronize_boundary_coordinates(), oomph::ThinLayerBrickOnTetMesh< ELEMENT >::ThinLayerBrickOnTetMesh(), and oomph::TwoDAnnularMesh< ELEMENT >::wrap_into_annular_shape().

◆ set_coordinates_on_boundary() [2/2]

virtual void oomph::Node::set_coordinates_on_boundary ( const unsigned &  b,
const Vector< double > &  boundary_zeta 
)
inlinevirtual

Set the vector of coordinates on mesh boundary b Broken virtual interface provides run-time error checking.

Definition at line 1428 of file nodes.h.

References set_coordinates_on_boundary().

◆ set_hanging_pt()

void oomph::Node::set_hanging_pt ( HangInfo *const &  hang_pt,
const int &  i 
)

Set the hanging data for the i-th value. (hang_pt=0 to make non-hanging)

Set the hanging data for the i-th value. If node is already hanging, simply overwrite the appropriate entry. If the node isn't hanging (because it might not be hanging geometrically), create the Vector of hanging pointers and make the other entries point to the node's geometric hanging data. Set hang_pt=0 to make entry explicitly non-hanging. Use Node::set_nonhanging() to unhang everything and clear up storage.

Definition at line 2068 of file nodes.cc.

References oomph::Data::constrain(), constrain_positions(), Hanging_pt, i, oomph::Data::nvalue(), oomph::Data::unconstrain(), and unconstrain_positions().

Referenced by oomph::TreeBasedRefineableMeshBase::complete_hanging_nodes(), oomph::RefineableQElement< 3 >::oc_hang_helper(), oomph::RefineableQElement< 2 >::quad_hang_helper(), oomph::Multi_domain_functions::recursively_add_masters_of_external_halo_node_to_storage(), oomph::Missing_masters_functions::recursively_add_masters_of_external_halo_node_to_storage(), and oomph::TreeBasedRefineableMeshBase::synchronise_hanging_nodes().

◆ set_non_obsolete()

void oomph::Node::set_non_obsolete ( )
inline

◆ set_nonhanging()

void oomph::Node::set_nonhanging ( )

◆ set_obsolete()

void oomph::Node::set_obsolete ( )
inline

◆ set_position_time_stepper()

void oomph::Node::set_position_time_stepper ( TimeStepper *const &  position_time_stepper_pt,
const bool &  preserve_existing_data 
)
virtual

Set a new position timestepper be resizing the appropriate storage.

Set a new position TimeStepper be resizing the appropriate storage. The current (zero) values will be unaffected, but all other entries will be set to zero.

Reimplemented in oomph::SolidNode.

Definition at line 1752 of file nodes.cc.

References ndim(), nposition_type(), oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, position_time_stepper_pt(), t, and X_position.

◆ unconstrain_positions()

virtual void oomph::Node::unconstrain_positions ( )
inlinevirtual

Unconstrain the positions when the node is made non-hanging Empty virtual function that is overloaded in SolidNodes.

Reimplemented in oomph::SolidNode.

Definition at line 1349 of file nodes.h.

Referenced by set_hanging_pt(), and set_nonhanging().

◆ unpin_all()

virtual void oomph::Node::unpin_all ( )
inlinevirtual

The unpin_all() function must be overloaded by SolidNode, so we put the virtual interface here to avoid virtual functions in Data.

Reimplemented in oomph::SolidNode.

Definition at line 1204 of file nodes.h.

References oomph::Data::unpin_all().

Referenced by oomph::SolidNode::unpin_all().

◆ value() [1/5]

Vector<double> oomph::Node::value ( ) const
inline

Return vector of values calculated using value(vector).

Definition at line 1502 of file nodes.h.

References oomph::Data::nvalue().

Referenced by value().

◆ value() [2/5]

double oomph::Node::value ( const unsigned &  i) const

Return i-th value (dofs or pinned) at this node either directly or via hanging node representation. Note that this REDFINES the interface in Data Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called.

Return i-th value (free or pinned) at this node either directly or via hanging node representation.

Definition at line 2408 of file nodes.cc.

References hanging_pt(), i, is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), and raw_value().

Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::BoundaryNode< NODE_TYPE >::clear_copied_pointers(), oomph::AxisymmetricVolumeConstraintBoundingElement::contribution_to_volume_flux(), oomph::LinearisedQCrouzeixRaviartElement::copy_efunction_to_normalisation(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::copy_onto_original_mesh(), oomph::Mesh::delete_all_external_storage(), oomph::ClampedSlidingHermiteBeamBoundaryConditionElement::fill_in_contribution_to_residuals(), oomph::ClampedHermiteShellBoundaryConditionElement::fill_in_contribution_to_residuals(), oomph::HeightControlElement::fill_in_contribution_to_residuals(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::FSIImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(), oomph::ImposeImpenetrabilityElement< ELEMENT >::fill_in_generic_contribution_to_residuals_parall_lagr_multiplier(), oomph::ImposeParallelOutflowElement< ELEMENT >::fill_in_generic_contribution_to_residuals_parall_lagr_multiplier(), oomph::LinearisedAxisymPoroelasticBJS_FSIElement< FLUID_BULK_ELEMENT, POROELASTICITY_BULK_ELEMENT >::fill_in_generic_residual_contribution_axisym_poroelastic_fsi(), oomph::BiharmonicFluidBoundaryElement::fill_in_generic_residual_contribution_biharmonic_boundary(), oomph::LinearisedFSIAxisymmetricNStNoSlipBCElementElement< FLUID_BULK_ELEMENT, SOLID_BULK_ELEMENT >::fill_in_generic_residual_contribution_fsi_no_slip_axisym(), oomph::Problem::get_data_to_be_sent_during_load_balancing(), oomph::Problem::get_dofs(), oomph::Missing_masters_functions::get_required_master_nodal_information_helper(), oomph::Multi_domain_functions::get_required_master_nodal_information_helper(), oomph::Missing_masters_functions::get_required_nodal_information_helper(), oomph::Multi_domain_functions::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_load_balance_helper(), oomph::BiharmonicEquations< DIM >::interpolated_dudx(), oomph::DGFaceElement::interpolated_u(), oomph::FiniteElement::nodal_value(), oomph::LinearElasticitySmoothMesh< LINEAR_ELASTICITY_ELEMENT >::operator()(), oomph::PoissonSmoothMesh< POISSON_ELEMENT >::operator()(), oomph::ClampedHermiteShellBoundaryConditionElement::output(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::output(), oomph::FSIImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::output(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::RefineableImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::RefineableFSIImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::refineable_fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::square_of_l2_norm_of_error(), oomph::TimeStepper::time_derivative(), oomph::FluidInterfaceElement::u(), and oomph::BiharmonicElement< DIM >::u().

◆ value() [3/5]

double oomph::Node::value ( const unsigned &  t,
const unsigned &  i 
) const

Return i-th value at time level t (t=0: present, t>0: previous) either directly or via hanging node representation. Note that this REDEFINES the interface in Data Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called.

Return i-th value (free or pinned) at this node at time level t either directly or via hanging node representation.

Definition at line 2437 of file nodes.cc.

References hanging_pt(), i, is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), raw_value(), and t.

◆ value() [4/5]

void oomph::Node::value ( const unsigned &  t,
Vector< double > &  values 
) const

Compute Vector of values (dofs or pinned) in this data at time level t (t=0: present; t>0: previous). This interface explicitly takes the hanging status into account. Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called.

Compute Vector of values (dofs or pinned) at this node at time level t (t=0: present; t>0: previous) either directly or via hanging node representation.

Definition at line 2482 of file nodes.cc.

References i, oomph::Data::nvalue(), t, and value().

◆ value() [5/5]

void oomph::Node::value ( Vector< double > &  values) const

Compute Vector of values for the Data value taking the hanging node status into account. Note that this REDEFINES the interface in Data Thus, the present function will be called provided that it is accessed through a pointer to a node i.e. Node* node_pt->value() will take hanging information into account. If a pointer to a Node has been explicitly down-cast to a pointer to Data then the "wrong" (Data) version of the function will be called.

Compute Vector of values (dofs or pinned) at this Data object either directly or via hanging node representation.

Definition at line 2466 of file nodes.cc.

References i, oomph::Data::nvalue(), and value().

◆ x() [1/4]

double& oomph::Node::x ( const unsigned &  i)
inline

Return the i-th nodal coordinate.

Definition at line 1060 of file nodes.h.

References i, Nposition_type, x_gen_range_check(), and X_position.

Referenced by oomph::NonLinearElasticitySmoothMesh< ELEMENT >::actions_before_newton_solve(), oomph::RefineableTetgenMesh< ELEMENT >::adapt(), oomph::RefineableTriangleMesh< ELEMENT >::adapt(), oomph::Problem::adapt(), oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::RefineableTriangleMesh< ELEMENT >::add_non_delete_vertices_from_boundary_helper(), oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::algebraic_node_update(), oomph::TetMeshBase::assess_mesh_quality(), oomph::Newmark< NSTEPS >::assign_initial_data_values(), oomph::GeneralisedElement::assign_local_eqn_numbers(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::backup(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableSolidQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::RefineableSolidQElement< 2 >::build(), oomph::GmshTetMesh< ELEMENT >::build_from_scaffold(), oomph::GeompackQuadMesh< ELEMENT >::build_from_scaffold(), oomph::SimpleCubicTetMesh< ELEMENT >::build_from_scaffold(), oomph::TetgenMesh< ELEMENT >::build_from_scaffold(), oomph::TriangleMesh< ELEMENT >::build_from_scaffold(), oomph::QuadFromTriangleMesh< ELEMENT >::build_from_scaffold(), oomph::ExtrudedCubeMeshFromQuadMesh< ELEMENT >::build_mesh(), oomph::FishMesh< ELEMENT >::build_mesh(), oomph::BrickFromTetMesh< ELEMENT >::build_mesh(), oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::BDF< NSTEPS >::calculate_predicted_positions(), oomph::TR::calculate_predicted_positions(), oomph::Mesh::check_for_repeated_nodes(), oomph::Mesh::check_halo_schemes(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::check_integrity(), oomph::RefineableQElement< 1 >::check_integrity(), oomph::RefineableQElement< 2 >::check_integrity(), oomph::FiniteElement::check_jacobian(), oomph::RefineableTriangleMesh< ELEMENT >::compute_area_target(), oomph::TriangleMesh< ELEMENT >::compute_boundary_segments_connectivity_and_initial_zeta_values(), oomph::RefineableTriangleMesh< ELEMENT >::compute_global_node_names_and_shared_nodes(), oomph::TriangleMesh< ELEMENT >::compute_holes_left_by_halo_elements_helper(), oomph::RefineableTriangleMesh< ELEMENT >::compute_shared_node_degree_helper(), oomph::RefineableTetMeshBase::compute_volume_target(), oomph::Multi_domain_functions::construct_new_external_halo_master_node_helper(), oomph::Missing_masters_functions::construct_new_external_halo_master_node_helper(), oomph::Multi_domain_functions::construct_new_external_halo_node_helper(), oomph::Missing_masters_functions::construct_new_external_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_node_load_balance_helper(), oomph::RefineableTriangleMesh< ELEMENT >::create_adjacency_matrix_new_shared_edges_helper(), oomph::TwoDimensionalPMLHelper::create_bottom_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_right_pml_mesh(), oomph::GmshTetScaffoldMesh::create_mesh_from_msh_file(), oomph::RefineableTriangleMesh< ELEMENT >::create_new_shared_boundaries(), oomph::RefineableTriangleMesh< ELEMENT >::create_sorted_face_mesh_representation(), oomph::TwoDimensionalPMLHelper::create_top_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_right_pml_mesh(), oomph::Mesh::delete_all_external_storage(), oomph::Mesh::distribute(), oomph::FSI_functions::doc_fsi(), oomph::Mesh::doc_mesh_distribution(), oomph::Mesh::doc_shared_nodes(), oomph::DummyErrorEstimator::DummyErrorEstimator(), dx_dt(), oomph::EighthSphereMesh< ELEMENT >::EighthSphereMesh(), oomph::FourierDecomposedHelmholtzDtNBoundaryElement< ELEMENT >::fill_in_generic_residual_contribution_fourier_decomposed_helmholtz_DtN_bc(), oomph::FSIDrivenCavityMesh< ELEMENT >::FSIDrivenCavityMesh(), oomph::FullCircleMesh< ELEMENT >::FullCircleMesh(), oomph::Problem::get_data_to_be_sent_during_load_balancing(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::RefineableElement::get_dresidual_dnodal_coordinates(), oomph::NavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::SpaceTimeNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableSpaceTimeNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::Z2ErrorEstimator::get_element_errors(), oomph::DummyErrorEstimator::get_element_errors(), oomph::RefineableTriangleMesh< ELEMENT >::get_face_mesh_representation(), oomph::Missing_masters_functions::get_required_master_nodal_information_helper(), oomph::Multi_domain_functions::get_required_master_nodal_information_helper(), oomph::Missing_masters_functions::get_required_nodal_information_helper(), oomph::Multi_domain_functions::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_helper(), oomph::RefineableTriangleMesh< ELEMENT >::get_required_nodal_information_load_balance_helper(), hanging_pt(), oomph::TriangleMesh< ELEMENT >::identify_boundary_segments_and_assign_initial_zeta_values(), oomph::HelmholtzMGPreconditioner< DIM >::maximum_edge_width(), oomph::NodeOrdering::node_global_position_comparison(), oomph::AlgebraicMesh::node_update(), oomph::Mesh::node_update(), oomph::MacroElementNodeUpdateNode::node_update(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_central_region(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::node_update_I(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_I(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::node_update_II(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_II(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::node_update_III(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_III(), oomph::AlgebraicFishMesh< ELEMENT >::node_update_in_body(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_central_box(), oomph::AlgebraicFishMesh< ELEMENT >::node_update_in_fin(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_lower_right_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_upper_left_box(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::node_update_IV(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_IV(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_IX(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_lower_right_region(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_upper_left_region(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_V(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_VI(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_VII(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::node_update_VIII(), oomph::RefineableQElement< 3 >::oc_hang_helper(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::oc_hang_helper(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::operator()(), oomph::LinearElasticitySmoothMesh< LINEAR_ELASTICITY_ELEMENT >::operator()(), oomph::PoissonSmoothMesh< POISSON_ELEMENT >::operator()(), oomph::DummyFaceElement< ELEMENT >::output(), output(), oomph::TriangleMesh< ELEMENT >::output_boundary_coordinates(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::output_integration_points(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::PRefineableQElement< 1, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), oomph::Tree::p_refine_if_required(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), oomph::MGSolver< DIM >::plot(), oomph::PMLQuadMeshBase< ELEMENT >::pml_locate_zeta(), position(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::position(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project(), oomph::RefineableQElement< 2 >::quad_hang_helper(), oomph::QuarterPipeMesh< ELEMENT >::QuarterPipeMesh(), oomph::QuarterTubeMesh< ELEMENT >::QuarterTubeMesh(), oomph::FiniteElement::raw_nodal_position(), oomph::TriangleMesh< ELEMENT >::re_assign_initial_zeta_values_for_internal_boundary(), oomph::TriangleMesh< ELEMENT >::re_scale_re_assigned_initial_zeta_values_for_internal_boundary(), oomph::Data::read(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::TreeBasedRefineableMeshBase::refine_as_in_reference_mesh(), oomph::Problem::remove_duplicate_data(), oomph::NonLinearElasticitySmoothMesh< ELEMENT >::reset(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::scale_basis(), oomph::TRaviartThomasDarcyElement< ORDER >::scale_basis(), oomph::TPoroelasticityElement< ORDER >::scale_basis(), oomph::Mesh::scale_mesh(), oomph::TriangleMesh< ELEMENT >::select_boundary_face_elements(), oomph::AlgebraicNode::self_test(), oomph::RefineableTriangleMesh< ELEMENT >::send_boundary_node_info_of_shared_nodes(), oomph::MGSolver< DIM >::set_self_test_vector(), oomph::AlgebraicChannelWithLeafletMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicCylinderWithFlagMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::setup_algebraic_node_update(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::setup_algebraic_node_update(), oomph::TetMeshBase::setup_boundary_coordinates(), oomph::XdaTetMesh< ELEMENT >::setup_boundary_coordinates(), oomph::UnstructuredTwoDMeshGeometryBase::setup_boundary_coordinates(), oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::setup_gamma(), oomph::HelmholtzDtNMesh< ELEMENT >::setup_gamma(), oomph::SimpleCubicScaffoldTetMesh::SimpleCubicScaffoldTetMesh(), oomph::RefineableTetgenMesh< ELEMENT >::snap_nodes_onto_boundary(), oomph::RefineableTriangleMesh< ELEMENT >::snap_nodes_onto_boundary(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), oomph::UnstructuredTwoDMeshGeometryBase::snap_nodes_onto_geometric_objects(), oomph::TetMeshBase::snap_to_quadratic_surface(), oomph::RefineableTriangleMesh< ELEMENT >::sort_nodes_on_shared_boundaries(), oomph::TwoDimensionalPMLHelper::sorter_bottom_boundary(), oomph::TwoDimensionalPMLHelper::sorter_left_boundary(), oomph::TwoDimensionalPMLHelper::sorter_right_boundary(), oomph::TwoDimensionalPMLHelper::sorter_top_boundary(), oomph::ChannelSpineMesh< ELEMENT >::spine_node_update(), oomph::HorizontalSingleLayerSpineMesh< ELEMENT >::spine_node_update(), oomph::SingleLayerCubicSpineMesh< ELEMENT >::spine_node_update(), oomph::SingleLayerSpineMesh< ELEMENT >::spine_node_update(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_channel(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_film_lower(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_film_upper(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_horizontal_transition_lower(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_horizontal_transition_upper(), oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update_lower(), oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update_upper(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_vertical_transition_lower(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_vertical_transition_upper(), oomph::TetMeshBase::split_elements_in_corners(), oomph::Tree::split_if_required(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::strain_rate(), oomph::TreeBasedRefineableMeshBase::synchronise_nonhanging_nodes(), oomph::TriangleMesh< ELEMENT >::synchronize_boundary_coordinates(), oomph::BDF< NSTEPS >::temporal_error_in_position(), oomph::TR::temporal_error_in_position(), oomph::TetMeshVertex::TetMeshVertex(), oomph::ThinLayerBrickOnTetMesh< ELEMENT >::ThinLayerBrickOnTetMesh(), oomph::TubeMesh< ELEMENT >::TubeMesh(), oomph::RefineableTriangleMesh< ELEMENT >::update_open_curve_after_restart(), oomph::RefineableTriangleMesh< ELEMENT >::update_open_curve_using_elements_area(), oomph::RefineableTriangleMesh< ELEMENT >::update_open_curve_using_face_mesh(), oomph::RefineableTriangleMesh< ELEMENT >::update_polygon_after_restart(), oomph::RefineableTriangleMesh< ELEMENT >::update_polygon_using_elements_area(), oomph::RefineableTriangleMesh< ELEMENT >::update_polygon_using_face_mesh(), oomph::TriangleMesh< ELEMENT >::update_triangulateio(), oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh(), and oomph::TwoDAnnularMesh< ELEMENT >::wrap_into_annular_shape().

◆ x() [2/4]

const double& oomph::Node::x ( const unsigned &  i) const
inline

Return the i-th nodal coordinate (const version).

Definition at line 1069 of file nodes.h.

References i, Nposition_type, x_gen_range_check(), and X_position.

◆ x() [3/4]

double& oomph::Node::x ( const unsigned &  t,
const unsigned &  i 
)
inline

Return the position x(i) at previous timestep t (t=0: present; t>0 previous timestep).

Definition at line 1079 of file nodes.h.

References i, Nposition_type, t, x_gen_range_check(), and X_position.

◆ x() [4/4]

const double& oomph::Node::x ( const unsigned &  t,
const unsigned &  i 
) const
inline

Return the position x(i) at previous timestep t (t=0: present; t>0 previous timestep) (const version)

Definition at line 1089 of file nodes.h.

References i, Nposition_type, t, x_gen_range_check(), and X_position.

◆ x_gen() [1/4]

double& oomph::Node::x_gen ( const unsigned &  k,
const unsigned &  i 
)
inline

Reference to the generalised position x(k,i). ‘Type’: k; Coordinate direction: i.

Definition at line 1126 of file nodes.h.

References i, Nposition_type, x_gen_range_check(), and X_position.

Referenced by oomph::Steady< NSTEPS >::assign_initial_positions_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_positions_impulsive(), oomph::BDF< NSTEPS >::assign_initial_positions_impulsive(), oomph::RefineableSolidElement::assign_solid_hanging_local_eqn_numbers(), oomph::SolidFiniteElement::assign_solid_local_eqn_numbers(), dx_gen_dt(), oomph::ClampedHermiteShellBoundaryConditionElement::fill_in_contribution_to_residuals(), oomph::BiharmonicFluidBoundaryElement::fill_in_generic_residual_contribution_biharmonic_boundary(), oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineablePseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd(), position_gen(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project(), oomph::FiniteElement::raw_nodal_position_gen(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::restore_positions(), oomph::SolidMesh::set_lagrangian_nodal_coordinates(), oomph::Problem::set_pinned_values_to_zero(), oomph::HermiteQuadMesh< ELEMENT >::set_position_of_node(), oomph::IMRBase::shift_time_positions(), oomph::Steady< NSTEPS >::shift_time_positions(), oomph::Newmark< NSTEPS >::shift_time_positions(), oomph::NewmarkBDF< NSTEPS >::shift_time_positions(), oomph::BDF< NSTEPS >::shift_time_positions(), and oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::store_positions().

◆ x_gen() [2/4]

const double& oomph::Node::x_gen ( const unsigned &  k,
const unsigned &  i 
) const
inline

Reference to the generalised position x(k,i). ‘Type’: k; Coordinate direction: i (const version).

Definition at line 1136 of file nodes.h.

References i, Nposition_type, x_gen_range_check(), and X_position.

◆ x_gen() [3/4]

double& oomph::Node::x_gen ( const unsigned &  t,
const unsigned &  k,
const unsigned &  i 
)
inline

Reference to the generalised position x(k,i) at the previous timestep [t=0: present]. ‘Type’: k; Coordinate direction: i.

Definition at line 1146 of file nodes.h.

References i, Nposition_type, t, x_gen_range_check(), and X_position.

◆ x_gen() [4/4]

const double& oomph::Node::x_gen ( const unsigned &  t,
const unsigned &  k,
const unsigned &  i 
) const
inline

Reference to the generalised position x(k,i) at the previous timestep [t=0: present]. ‘Type’: k; Coordinate direction: i. (const version)

Definition at line 1157 of file nodes.h.

References i, Nposition_type, t, x_gen_range_check(), and X_position.

◆ x_gen_range_check()

void oomph::Node::x_gen_range_check ( const unsigned &  t,
const unsigned &  k,
const unsigned &  i 
) const
protected

Private function to check that the arguemnts to the position functions are in range.

///////////////////////////////////////////////////////////////////////

/////////////////////////////////////////////////////////////////////// Private function to check that the arguments are within the range of the stored coordinates, position types and time history values.

Definition at line 1528 of file nodes.cc.

References i, Ndim, Nposition_type, oomph::TimeStepper::ntstorage(), Position_time_stepper_pt, and t.

Referenced by x(), and x_gen().

◆ x_position_pt()

double* oomph::Node::x_position_pt ( const unsigned &  i)
inlineprotected

Direct access to the pointer to the i-th stored coordinate data.

Definition at line 958 of file nodes.h.

References i, and X_position.

◆ x_pt()

double* oomph::Node::x_pt ( const unsigned &  t,
const unsigned &  i 
)
inline

Direct access to the i-th coordinate at time level t (t=0: present; t>0: previous)

Definition at line 1181 of file nodes.h.

References i, Nposition_type, t, and X_position.

Friends And Related Function Documentation

◆ BoundaryNodeBase

friend class BoundaryNodeBase
friend

to construct periodic Nodes

Definition at line 914 of file nodes.h.

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const Node d 
)
friend

Output operator: output location and all values at all times, along with any extra information stored for the timestepper.

Definition at line 373 of file nodes.cc.

Member Data Documentation

◆ Aux_node_update_fct_pt

AuxNodeUpdateFctPt oomph::Node::Aux_node_update_fct_pt
protected

Pointer to auxiliary update function – this can be used to update any nodal values following the update of the nodal position. This is needed e.g. to update the no-slip condition on moving boundaries.

Definition at line 967 of file nodes.h.

Referenced by has_auxiliary_node_update_fct_pt(), oomph::AlgebraicNode::node_update(), oomph::MacroElementNodeUpdateNode::node_update(), oomph::SpineNode::node_update(), perform_auxiliary_node_update_fct(), and set_auxiliary_node_update_fct_pt().

◆ Hanging_pt

HangInfo** oomph::Node::Hanging_pt
protected

C-style array of pointers to hanging node info. It's set to NULL if the node isn't hanging. The first entry (0) is the geometric hanging node data. The remaining entries correspond to the hanging data for the other values stored at the node. Usually, these entries will be the same as the geometric hanging node data represented by Hanging_pt[0], but this is not necessarily the case; e.g. the pressure in Taylor Hood has different hanging node data from the velocities.

Definition at line 942 of file nodes.h.

Referenced by hanging_pt(), is_hanging(), resize(), set_hanging_pt(), set_nonhanging(), and ~Node().

◆ Ndim

unsigned oomph::Node::Ndim
protected

Eulerian dimension of the node.

Definition at line 945 of file nodes.h.

Referenced by copy(), dump(), ndim(), read(), oomph::SolidNode::SolidNode(), and x_gen_range_check().

◆ No_independent_position

unsigned oomph::Node::No_independent_position = 10
static

Static "Magic number" used to indicate that there is no independent position in a periodic node.

Static "Magic number" passed as independent_position when there is no independent position in the periodic node. For example, in a periodic mesh.

Definition at line 972 of file nodes.h.

◆ Nposition_type

unsigned oomph::Node::Nposition_type
protected

Number of coordinate types used in the mapping between local and global coordinates (e.g. 1 for Lagrange-type elements; 2 for 1D Hermite elements; 4 for 2D Hermite elements, etc).

Definition at line 951 of file nodes.h.

Referenced by copy(), dump(), nposition_type(), oomph::SolidNode::pin_position(), oomph::SolidNode::position_eqn_number(), oomph::SolidNode::position_is_a_copy(), oomph::SolidNode::position_is_pinned(), read(), oomph::SolidNode::SolidNode(), oomph::SolidNode::unpin_position(), x(), x_gen(), x_gen_range_check(), and x_pt().

◆ Obsolete

bool oomph::Node::Obsolete
protected

Flag to indicate that the Node has become obsolete — usually during mesh refinement process.

Definition at line 955 of file nodes.h.

Referenced by is_obsolete(), set_non_obsolete(), and set_obsolete().

◆ Position_time_stepper_pt

TimeStepper* oomph::Node::Position_time_stepper_pt
protected

◆ X_position

double** oomph::Node::X_position
protected

Array of pointers to the data holding the Eulerian positions. The storage format must be the same as the internal data storage so that we can implement the functions x() in generality here without the need for virtual functions. The first index will be a flat array of position types and coordinates and the second will be the number of time history values at each position type.

Definition at line 929 of file nodes.h.

Referenced by add_values_to_vector(), copy(), dump(), Node(), read(), read_values_from_vector(), oomph::SolidNode::set_external_variable_position_pt(), set_position_time_stepper(), oomph::SolidNode::set_position_time_stepper(), oomph::SolidNode::SolidNode(), x(), x_gen(), x_position_pt(), x_pt(), ~Node(), and oomph::SolidNode::~SolidNode().


The documentation for this class was generated from the following files: