Displacement control element: In the "normal" formulation of solid mechanics problems, the external load is given and the displacement throughout the solid body is computed. For highly nonlinear problems it is sometimes helpful to re-formulate the problem by prescribing the position of a selected control point and treating the (scalar) load level required to achieve this deformation as an unknown. As an example consider the buckling of pressure-loaded, thin-walled elastic shells. The load-displacement characteristics of such structures tend to be highly nonlinear and bifurcations from the structure's pre-buckling state often occur via sub-critical bifurcations. If we have some a-priori knowledge of the expected deformation (for example, during the non-axisymmetric buckling of a circular cylindrical shell certain material points will be displaced radially inwards), it is advantageous to prescribe the radial displacement of a carefully selected control point and treat the external pressure as an unknown. More...
#include <displacement_control_element.h>
Public Member Functions | |
DisplacementControlElement (SolidFiniteElement *controlled_element_pt, const Vector< double > &controlled_point, const unsigned &controlled_direction, double *control_position_value_pt, Data *displacement_control_load_pt) | |
Constructor. Pass: More... | |
DisplacementControlElement (SolidFiniteElement *controlled_element_pt, const Vector< double > &controlled_point, const unsigned &controlled_direction, double *control_position_value_pt) | |
Constructor. Pass: More... | |
DisplacementControlElement (const DisplacementControlElement &)=delete | |
Broken copy constructor. More... | |
void | operator= (const DisplacementControlElement &)=delete |
Broken assignment operator. More... | |
Data * | displacement_control_load_pt () const |
Pointer to Data object whose one-and-only value represents the load that is adjusted to allow displacement control. More... | |
void | assign_additional_local_eqn_numbers () |
Store local equation number of displacement control equation. More... | |
void | fill_in_contribution_to_residuals (Vector< double > &residuals) |
Add the element's contribution to its residual vector: The displacement constraint. [Note: Jacobian is computed automatically by finite-differencing]. More... | |
unsigned | ndof_types () const |
The number of "DOF" that degrees of freedom in this element are sub-divided into: Just the control pressure. More... | |
void | get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const |
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contains the global equation number of the unknown, while the second one contains the number of the "DOF type" that this unknown is associated with. (Function can obviously only be called if the equation numbering scheme has been set up.) The only dof this element is in charge of is the control load, provided it's been created as internal Data. More... | |
Public Member Functions inherited from oomph::GeneralisedElement | |
GeneralisedElement() | GeneralisedElement (const GeneralisedElement &)=delete |
Constructor: Initialise all pointers and all values to zero. More... | |
void | operator= (const GeneralisedElement &)=delete |
Broken assignment operator. More... | |
Data *& | internal_data_pt (const unsigned &i) |
Return a pointer to i-th internal data object. More... | |
Data *const & | internal_data_pt (const unsigned &i) const |
Return a pointer to i-th internal data object (const version) More... | |
Data *& | external_data_pt (const unsigned &i) |
Return a pointer to i-th external data object. More... | |
Data *const & | external_data_pt (const unsigned &i) const |
Return a pointer to i-th external data object (const version) More... | |
unsigned long | eqn_number (const unsigned &ieqn_local) const |
Return the global equation number corresponding to the ieqn_local-th local equation number. More... | |
int | local_eqn_number (const unsigned long &ieqn_global) const |
Return the local equation number corresponding to the ieqn_global-th global equation number. Returns minus one (-1) if there is no local degree of freedom corresponding to the chosen global equation number. More... | |
unsigned | add_external_data (Data *const &data_pt, const bool &fd=true) |
Add a (pointer to an) external data object to the element and return its index (i.e. the index required to obtain it from the access function external_data_pt(...) . The optional boolean flag indicates whether the data should be included in the general finite-difference loop when calculating the jacobian. The default value is true, i.e. the data will be included in the finite-differencing. More... | |
bool | external_data_fd (const unsigned &i) const |
Return the status of the boolean flag indicating whether the external data is included in the finite difference loop. More... | |
void | exclude_external_data_fd (const unsigned &i) |
Set the boolean flag to exclude the external datum from the the finite difference loop when computing the jacobian matrix. More... | |
void | include_external_data_fd (const unsigned &i) |
Set the boolean flag to include the external datum in the the finite difference loop when computing the jacobian matrix. More... | |
void | flush_external_data () |
Flush all external data. More... | |
void | flush_external_data (Data *const &data_pt) |
Flush the object addressed by data_pt from the external data array. More... | |
unsigned | ninternal_data () const |
Return the number of internal data objects. More... | |
unsigned | nexternal_data () const |
Return the number of external data objects. More... | |
unsigned | ndof () const |
Return the number of equations/dofs in the element. More... | |
void | dof_vector (const unsigned &t, Vector< double > &dof) |
Return the vector of dof values at time level t. More... | |
void | dof_pt_vector (Vector< double * > &dof_pt) |
Return the vector of pointers to dof values. More... | |
void | set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data) |
Set the timestepper associated with the i-th internal data object. More... | |
void | assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt) |
Assign the global equation numbers to the internal Data. The arguments are the current highest global equation number (which will be incremented) and a Vector of pointers to the global variables (to which any unpinned values in the internal Data are added). More... | |
void | describe_dofs (std::ostream &out, const std::string ¤t_string) const |
Function to describe the dofs of the element. 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 | describe_local_dofs (std::ostream &out, const std::string ¤t_string) const |
Function to describe the local dofs of the element. 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... | |
void | add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt) |
Add pointers to the internal data values to map indexed by the global equation number. More... | |
void | add_internal_data_values_to_vector (Vector< double > &vector_of_values) |
Add all internal data and time history values to the vector in the internal storage order. More... | |
void | read_internal_data_values_from_vector (const Vector< double > &vector_of_values, unsigned &index) |
Read all internal data and time history values 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... | |
void | add_internal_eqn_numbers_to_vector (Vector< long > &vector_of_eqn_numbers) |
Add all equation numbers associated with internal data to the vector in the internal storage order. More... | |
void | read_internal_eqn_numbers_from_vector (const Vector< long > &vector_of_eqn_numbers, unsigned &index) |
Read all equation numbers associated with internal data 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... | |
virtual void | assign_local_eqn_numbers (const bool &store_local_dof_pt) |
Setup the arrays of local equation numbers for the element. If the optional boolean argument is true, then pointers to the associated degrees of freedom are stored locally in the array Dof_pt. More... | |
virtual void | complete_setup_of_dependencies () |
Complete the setup of any additional dependencies that the element may have. Empty virtual function that may be overloaded for specific derived elements. Used, e.g., for elements with algebraic node update functions to determine the "geometric
Data", i.e. the Data that affects the element's shape. This function is called (for all elements) at the very beginning of the equation numbering procedure to ensure that all dependencies are accounted for. More... | |
virtual void | get_residuals (Vector< double > &residuals) |
Calculate the vector of residuals of the equations in the element. By default initialise the vector to zero and then call the fill_in_contribution_to_residuals() function. Note that this entire function can be overloaded if desired. More... | |
virtual void | get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Calculate the elemental Jacobian matrix "d equation / d
variable". More... | |
virtual void | get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix) |
Calculate the residuals and the elemental "mass" matrix, the matrix that multiplies the time derivative terms in a problem. More... | |
virtual void | get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
Calculate the residuals and jacobian and elemental "mass" matrix, the matrix that multiplies the time derivative terms. More... | |
virtual void | get_dresiduals_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam) |
Calculate the derivatives of the residuals with respect to a parameter. More... | |
virtual void | get_djacobian_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
Calculate the derivatives of the elemental Jacobian matrix and residuals with respect to a parameter. More... | |
virtual void | get_djacobian_and_dmass_matrix_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam) |
Calculate the derivatives of the elemental Jacobian matrix mass matrix and residuals with respect to a parameter. More... | |
virtual void | get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product) |
Calculate the product of the Hessian (derivative of Jacobian with respect to all variables) an eigenvector, Y, and other specified vectors, C (d(J_{ij})/d u_{k}) Y_{j} C_{k}. More... | |
virtual void | get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product) |
Return the vector of inner product of the given pairs of history values. More... | |
virtual void | get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector) |
Compute the vectors that when taken as a dot product with other history values give the inner product over the element. More... | |
virtual unsigned | self_test () |
Self-test: Have all internal values been classified as pinned/unpinned? Return 0 if OK. More... | |
virtual void | compute_norm (Vector< double > &norm) |
Compute norm of solution – broken virtual can be overloaded by element writer to implement whatever norm is desired for the specific element. More... | |
virtual void | compute_norm (double &norm) |
Compute norm of solution – broken virtual can be overloaded by element writer to implement whatever norm is desired for the specific element. More... | |
void | set_halo (const unsigned &non_halo_proc_ID) |
Label the element as halo and specify processor that holds non-halo counterpart. More... | |
void | set_nonhalo () |
Label the element as not being a halo. More... | |
bool | is_halo () const |
Is this element a halo? More... | |
int | non_halo_proc_ID () |
ID of processor ID that holds non-halo counterpart of halo element; negative if not a halo. More... | |
void | set_must_be_kept_as_halo () |
Insist that this element be kept as a halo element during a distribute? More... | |
void | unset_must_be_kept_as_halo () |
Do not insist that this element be kept as a halo element during distribution. More... | |
bool | must_be_kept_as_halo () const |
Test whether the element must be kept as a halo element. More... | |
Protected Attributes | |
Data * | Displacement_control_load_pt |
Pointer to Data item whose one-and-only value contains the load value that is being adjusted to allow displacement control. More... | |
double * | Control_position_value_pt |
Pointer to the value that stores the prescribed coordinate of the control point. More... | |
unsigned | Controlled_direction |
Coordinate direction in which the displacement of the control point is controlled. More... | |
SolidFiniteElement * | Controlled_element_pt |
Pointer to SolidFiniteElement at which control displacement is applied. More... | |
Vector< double > | Controlled_point |
Vector of local coordinates of point at which control displacement is applied. More... | |
bool | Load_data_created_internally |
Flag to indicate if load data was created internally or externally (and is therefore stored in the element's internal or external Data) More... | |
unsigned | Load_data_index |
In which component (in the vector of the element's internal or external Data) is the load stored? More... | |
int | Displ_ctrl_local_eqn |
Local equation number of the control-displacement equation. More... | |
Protected Attributes inherited from oomph::GeneralisedElement | |
int | Non_halo_proc_ID |
Non-halo processor ID for Data; -1 if it's not a halo. More... | |
bool | Must_be_kept_as_halo |
Does this element need to be kept as a halo element during a distribute? More... | |
Additional Inherited Members | |
Static Public Attributes inherited from oomph::GeneralisedElement | |
static bool | Suppress_warning_about_repeated_internal_data |
Static boolean to suppress warnings about repeated internal data. Defaults to false. More... | |
static bool | Suppress_warning_about_repeated_external_data = true |
Static boolean to suppress warnings about repeated external data. Defaults to true. More... | |
static double | Default_fd_jacobian_step = 1.0e-8 |
Double used for the default finite difference step in elemental jacobian calculations. More... | |
Protected Member Functions inherited from oomph::GeneralisedElement | |
unsigned | add_internal_data (Data *const &data_pt, const bool &fd=true) |
Add a (pointer to an) internal data object to the element and return the index required to obtain it from the access function internal_data_pt() . The boolean indicates whether the datum should be included in the general finite-difference loop when calculating the jacobian. The default value is true, i.e. the data will be included in the finite differencing. More... | |
bool | internal_data_fd (const unsigned &i) const |
Return the status of the boolean flag indicating whether the internal data is included in the finite difference loop. More... | |
void | exclude_internal_data_fd (const unsigned &i) |
Set the boolean flag to exclude the internal datum from the finite difference loop when computing the jacobian matrix. More... | |
void | include_internal_data_fd (const unsigned &i) |
Set the boolean flag to include the internal datum in the finite difference loop when computing the jacobian matrix. More... | |
void | clear_global_eqn_numbers () |
Clear the storage for the global equation numbers and pointers to dofs (if stored) More... | |
void | add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt) |
Add the contents of the queue global_eqn_numbers to the local storage for the local-to-global translation scheme. It is essential that the entries in the queue are added IN ORDER i.e. from the front. More... | |
virtual void | assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt) |
Assign the local equation numbers for the internal and external Data This must be called after the global equation numbers have all been assigned. It is virtual so that it can be overloaded by ElementWithExternalElements so that any external data from the external elements in included in the numbering scheme. If the boolean argument is true then pointers to the dofs will be stored in Dof_pt. More... | |
virtual void | assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt) |
Assign all the local equation numbering schemes that can be applied generically for the element. In most cases, this is the function that will be overloaded by inherited classes. It is required to ensure that assign_additional_local_eqn_numbers() can always be called after ALL other local equation numbering has been performed. The default for the GeneralisedElement is simply to call internal and external local equation numbering. If the boolean argument is true then pointers to the dofs will be stored in Dof_pt. More... | |
int | internal_local_eqn (const unsigned &i, const unsigned &j) const |
Return the local equation number corresponding to the j-th value stored at the i-th internal data. More... | |
int | external_local_eqn (const unsigned &i, const unsigned &j) |
Return the local equation number corresponding to the j-th value stored at the i-th external data. More... | |
void | fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
Calculate the contributions to the jacobian from the internal degrees of freedom using finite differences. This version of the function assumes that the residuals vector has already been calculated. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data. More... | |
void | fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
Calculate the contributions to the jacobian from the internal degrees of freedom using finite differences. This version computes the residuals vector before calculating the jacobian terms. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data. More... | |
void | fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
Calculate the contributions to the jacobian from the external degrees of freedom using finite differences. This version of the function assumes that the residuals vector has already been calculated. If the boolean argument is true, the finite differencing will be performed for all external data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data. More... | |
void | fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
Calculate the contributions to the jacobian from the external degrees of freedom using finite differences. This version computes the residuals vector before calculating the jacobian terms. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data. More... | |
virtual void | update_before_internal_fd () |
Function that is called before the finite differencing of any internal data. This may be overloaded to update any dependent data before finite differencing takes place. More... | |
virtual void | reset_after_internal_fd () |
Function that is call after the finite differencing of the internal data. This may be overloaded to reset any dependent variables that may have changed during the finite differencing. More... | |
virtual void | update_in_internal_fd (const unsigned &i) |
Function called within the finite difference loop for internal data after a change in any values in the i-th internal data object. More... | |
virtual void | reset_in_internal_fd (const unsigned &i) |
Function called within the finite difference loop for internal data after the values in the i-th external data object are reset. The default behaviour is to call the update function. More... | |
virtual void | update_before_external_fd () |
Function that is called before the finite differencing of any external data. This may be overloaded to update any dependent data before finite differencing takes place. More... | |
virtual void | reset_after_external_fd () |
Function that is call after the finite differencing of the external data. This may be overloaded to reset any dependent variables that may have changed during the finite differencing. More... | |
virtual void | update_in_external_fd (const unsigned &i) |
Function called within the finite difference loop for external data after a change in any values in the i-th external data object. More... | |
virtual void | reset_in_external_fd (const unsigned &i) |
Function called within the finite difference loop for external data after the values in the i-th external data object are reset. The default behaviour is to call the update function. More... | |
virtual void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this function will NOT initialise the residuals vector or the jacobian matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is to use finite differences to calculate the jacobian. More... | |
virtual void | fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix) |
Add the elemental contribution to the mass matrix matrix. and the residuals vector. Note that this function should NOT initialise the residuals vector or the mass matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is deliberately broken. More... | |
virtual void | fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
Add the elemental contribution to the jacobian matrix, mass matrix and the residuals vector. Note that this function should NOT initialise any entries. It must be called after the residuals vector and matrices have been initialised to zero. More... | |
virtual void | fill_in_contribution_to_dresiduals_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam) |
Add the elemental contribution to the derivatives of the residuals with respect to a parameter. This function should NOT initialise any entries and must be called after the entries have been initialised to zero The default implementation is to use finite differences to calculate the derivatives. More... | |
virtual void | fill_in_contribution_to_djacobian_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
Add the elemental contribution to the derivatives of the elemental Jacobian matrix and residuals with respect to a parameter. This function should NOT initialise any entries and must be called after the entries have been initialised to zero The default implementation is to use finite differences to calculate the derivatives. More... | |
virtual void | fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam) |
Add the elemental contribution to the derivative of the jacobian matrix, mass matrix and the residuals vector with respect to the passed parameter. Note that this function should NOT initialise any entries. It must be called after the residuals vector and matrices have been initialised to zero. More... | |
virtual void | fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product) |
Fill in contribution to the product of the Hessian (derivative of Jacobian with respect to all variables) an eigenvector, Y, and other specified vectors, C (d(J_{ij})/d u_{k}) Y_{j} C_{k}. More... | |
virtual void | fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product) |
Fill in the contribution to the inner products between given pairs of history values. More... | |
virtual void | fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector) |
Fill in the contributions to the vectors that when taken as dot product with other history values give the inner product over the element. More... | |
Static Protected Attributes inherited from oomph::GeneralisedElement | |
static DenseMatrix< double > | Dummy_matrix |
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case when only the residuals are being assembled. More... | |
static std::deque< double * > | Dof_pt_deque |
Static storage for deque used to add_global_equation_numbers when pointers to the dofs in each element are not required. More... | |
Displacement control element: In the "normal" formulation of solid mechanics problems, the external load is given and the displacement throughout the solid body is computed. For highly nonlinear problems it is sometimes helpful to re-formulate the problem by prescribing the position of a selected control point and treating the (scalar) load level required to achieve this deformation as an unknown. As an example consider the buckling of pressure-loaded, thin-walled elastic shells. The load-displacement characteristics of such structures tend to be highly nonlinear and bifurcations from the structure's pre-buckling state often occur via sub-critical bifurcations. If we have some a-priori knowledge of the expected deformation (for example, during the non-axisymmetric buckling of a circular cylindrical shell certain material points will be displaced radially inwards), it is advantageous to prescribe the radial displacement of a carefully selected control point and treat the external pressure as an unknown.
DisplacementControlElements
facilitate the use of such methods. They require the specification of
controlled_element_pt
, to a SolidFiniteElement
andcontrolled_point
which contains the local coordinates of the control point in that SolidFiniteElement
.controlled_direction
. in which the displacement is controlled.control_position_value_pt
, that specifies the desired value of the prescribed coordinate after the deformation (i.e. if controlled_direction=1
then *control_position_value_pt
specifies the coordinate of the control point in the deformed configuration.)The DisplacementControlElement
has two constructors:
Data
object whose one-and-only value contains the scalar load level that is "traded" for the displacement constraint. This is appropriate if the load Data
has already been created (and is included in the overall equation numbering procedure) by some other element. In that case the DisplacementControlElement
treats the (already existing) Data
object external Data
.Data
object (with a single value) is created by the constructor of the DisplacementControlElement
and stored in its internal Data
. Once the DisplacementControlElement
has been included in one of the Problem's
meshes, it is therefore automatically included in the equation numbering procedure. The (pointer to) the newly created Data
is accessible via the access function displacement_control_load_pt()
. It can be used to make make the unknown load level accessible to the load function that drives the deformation.Note: The element inherits from the BlockPreconditionableElementBase and can be used in the block-preconditioning context. The element is "in charge" of the control load (if it's been created internally) and classifies it as its one-and-only "DOF type"
Definition at line 105 of file displacement_control_element.h.
|
inline |
Constructor. Pass:
SolidFiniteElement
that contains the control pointThe load Data
is treated as external Data
for this element.
Definition at line 122 of file displacement_control_element.h.
References oomph::GeneralisedElement::add_external_data(), Controlled_element_pt, displacement_control_load_pt(), Displacement_control_load_pt, Load_data_created_internally, Load_data_index, oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::Data::nvalue(), and oomph::SolidNode::variable_position_pt().
|
inline |
Constructor. Pass:
SolidFiniteElement
that contains the control pointThe pointer to a Data item whose one-and-only value contains the load value that is being adjusted to allow displacement control is created internally (and stored in the element's internal Data
. It is accessible (for use the load function) via the access function displacement_control_load_pt()
Definition at line 177 of file displacement_control_element.h.
References oomph::GeneralisedElement::add_external_data(), oomph::GeneralisedElement::add_internal_data(), Controlled_element_pt, Displacement_control_load_pt, Load_data_created_internally, Load_data_index, oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), and oomph::SolidNode::variable_position_pt().
|
delete |
Broken copy constructor.
|
inlinevirtual |
Store local equation number of displacement control equation.
Reimplemented from oomph::GeneralisedElement.
Definition at line 226 of file displacement_control_element.h.
References Displ_ctrl_local_eqn, oomph::GeneralisedElement::external_local_eqn(), oomph::GeneralisedElement::internal_local_eqn(), Load_data_created_internally, and Load_data_index.
|
inline |
Pointer to Data object whose one-and-only value represents the load that is adjusted to allow displacement control.
Definition at line 219 of file displacement_control_element.h.
References Displacement_control_load_pt.
Referenced by DisplacementControlElement().
|
inlinevirtual |
Add the element's contribution to its residual vector: The displacement constraint. [Note: Jacobian is computed automatically by finite-differencing].
Reimplemented from oomph::GeneralisedElement.
Definition at line 248 of file displacement_control_element.h.
References Control_position_value_pt, Controlled_direction, Controlled_element_pt, Controlled_point, Displ_ctrl_local_eqn, and oomph::FiniteElement::interpolated_x().
|
inlinevirtual |
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contains the global equation number of the unknown, while the second one contains the number of the "DOF type" that this unknown is associated with. (Function can obviously only be called if the equation numbering scheme has been set up.) The only dof this element is in charge of is the control load, provided it's been created as internal Data.
Reimplemented from oomph::GeneralisedElement.
Definition at line 274 of file displacement_control_element.h.
References Displ_ctrl_local_eqn, oomph::GeneralisedElement::eqn_number(), Load_data_created_internally, and oomph::GeneralisedElement::local_eqn_number().
|
inlinevirtual |
The number of "DOF" that degrees of freedom in this element are sub-divided into: Just the control pressure.
Reimplemented from oomph::GeneralisedElement.
Definition at line 261 of file displacement_control_element.h.
|
delete |
Broken assignment operator.
|
protected |
Pointer to the value that stores the prescribed coordinate of the control point.
Definition at line 308 of file displacement_control_element.h.
Referenced by fill_in_contribution_to_residuals().
|
protected |
Coordinate direction in which the displacement of the control point is controlled.
Definition at line 312 of file displacement_control_element.h.
Referenced by fill_in_contribution_to_residuals().
|
protected |
Pointer to SolidFiniteElement at which control displacement is applied.
Definition at line 315 of file displacement_control_element.h.
Referenced by DisplacementControlElement(), and fill_in_contribution_to_residuals().
|
protected |
Vector of local coordinates of point at which control displacement is applied.
Definition at line 319 of file displacement_control_element.h.
Referenced by fill_in_contribution_to_residuals().
|
protected |
Local equation number of the control-displacement equation.
Definition at line 331 of file displacement_control_element.h.
Referenced by assign_additional_local_eqn_numbers(), fill_in_contribution_to_residuals(), and get_dof_numbers_for_unknowns().
|
protected |
Pointer to Data item whose one-and-only value contains the load value that is being adjusted to allow displacement control.
Definition at line 304 of file displacement_control_element.h.
Referenced by displacement_control_load_pt(), and DisplacementControlElement().
|
protected |
Flag to indicate if load data was created internally or externally (and is therefore stored in the element's internal or external Data)
Definition at line 324 of file displacement_control_element.h.
Referenced by assign_additional_local_eqn_numbers(), DisplacementControlElement(), and get_dof_numbers_for_unknowns().
|
protected |
In which component (in the vector of the element's internal or external Data) is the load stored?
Definition at line 328 of file displacement_control_element.h.
Referenced by assign_additional_local_eqn_numbers(), and DisplacementControlElement().