OGS
ProcessLib::ReleaseNodalForce Class Referencefinal

Detailed Description

Boundary condition for simulating excavation using the release nodal force approach.

This class implements a boundary condition that applies a time-dependent, released nodal force to nodes on an exposed surface for excavation simulations.

The initial state assumes a non-equilibrium stress \( \sigma_0 \). In the finite element method, the nodal force is given by:

\[ \mathbf{b} = \int \left( \text{B}^T (\mathbf{\sigma} - \mathbf{\sigma}_0) + (\mathbf{f} - \mathbf{f}_0) \right) N \, \mathrm{d}\Omega + \int_{\Gamma_q} (\boldsymbol{\sigma}-\boldsymbol{\sigma}_0)\cdot \mathbf n \mathrm{d}\Gamma \]

where:

  • \( \text{B} \) is the strain-displacement matrix,
  • \( \mathbf{\sigma} \) is the current total stress,
  • \( \mathbf{\sigma}_0 \) is the initial total stress,
  • \( \mathbf{f} \) is the current body force,
  • \( \mathbf{f}_0 \) is the initial body force,
  • \(\int_{\Gamma_q}\) is the boundary where the traction condition is applied, and \( \mathbf n\) is the outward normal of the boundary,
  • \( N \) is the shape function,
  • \( \Omega \) is the domain of integration.

After excavation, the stress and body force inside the excavated domain vanish, leaving non-zero nodal forces at the exposed surface nodes. These are computed as:

\[ \mathbf{b}_0 = -\int \left( \text{B}^T \mathbf{\sigma}_0 + \mathbf{f}_0 N \right) \, \mathrm{d}\Omega - \int_{\Gamma_q} \boldsymbol{\sigma}_0 \cdot \mathbf n \mathrm{d}\Gamma \]

where \(\Omega\) is the remaining domain.

The elements of \( \mathbf{b}_0 \) corresponding to the exposed surface nodes define the released nodal force vector:

\[ \mathbf{f}_\text{r} := (\mathbf{b}_0)_i, \quad i \in \text{exposed surface nodes}. \]

\(\Omega\) can be the excavated domain, which leads to the negative \( \mathbf{f}_\text{r}\).

To simulate excavations under an assumption of gradual release of these forces, the boundary condition applies the released nodal force vector to the global right-hand side (RHS) vector b, scaled by a time- and position-dependent release parameter \( g(t, \mathbf{x}) \):

\[ \mathbf{b} = \mathbf{b} + \mathbf{f}_\text{r} \cdot g(t, \mathbf{x}) \]

The release parameter should be a monotonically decreasing function, representing the progressive removal of support over time, e.g., \( g(0, \mathbf{x}) = 1 \) and \( g(t_e, \mathbf{x}) = 0 \), where \( t_e \) is the end time of excavation, and \( \frac{\partial g}{\partial t} < 0 \).

This boundary condition is particularly useful for modeling staged excavations or similar processes where loads are released in a controlled manner over time.

Note
Setting compensate_non_equilibrium_initial_residuum to true in the process variable configuration is required when using this boundary condition to ensure that the initial non-equilibrium stress state is properly accounted for.

Definition at line 100 of file ReleaseNodalForce.h.

#include <ReleaseNodalForce.h>

Inheritance diagram for ProcessLib::ReleaseNodalForce:
[legend]
Collaboration diagram for ProcessLib::ReleaseNodalForce:
[legend]

Public Member Functions

 ReleaseNodalForce (int const variable_id, MeshLib::Mesh const &boundary_mesh, std::unique_ptr< NumLib::LocalToGlobalIndexMap > &dof_table, ParameterLib::Parameter< double > const &time_decay_parameter)
 Constructs a released nodal force boundary condition.
 
void set (GlobalVector const *r_neq)
 
void applyNaturalBC (const double t, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix *K, GlobalVector &b, GlobalMatrix *Jac) override
 Applies the released nodal force boundary condition. This method scales the nodal forces by the release parameter and applies them to the global RHS vector b.
 
- Public Member Functions inherited from ProcessLib::BoundaryCondition
virtual void getEssentialBCValues (const double, GlobalVector const &, NumLib::IndexValueVector< GlobalIndexType > &) const
 Writes the values of essential BCs to bc_values.
 
virtual void preTimestep (const double, std::vector< GlobalVector * > const &, int const)
 
virtual void postTimestep (const double, std::vector< GlobalVector * > const &, int const)
 
virtual ~BoundaryCondition ()=default
 

Private Attributes

int const variable_id_
 
MeshLib::Mesh const & boundary_mesh_
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap > const dof_table_
 
ParameterLib::Parameter< double > const & time_decay_parameter_
 
std::vector< double > initial_release_nodal_force_
 
std::vector< GlobalIndexTypeglobal_indices_
 
std::vector< MeshLib::Node const * > boundary_nodes_
 

Constructor & Destructor Documentation

◆ ReleaseNodalForce()

ProcessLib::ReleaseNodalForce::ReleaseNodalForce ( int const variable_id,
MeshLib::Mesh const & boundary_mesh,
std::unique_ptr< NumLib::LocalToGlobalIndexMap > & dof_table,
ParameterLib::Parameter< double > const & time_decay_parameter )
explicit

Constructs a released nodal force boundary condition.

Parameters
variable_idThe ID of the variable to which the nodal forces are applied.
boundary_meshThe mesh representing the boundary where the nodal forces are applied.
dof_tableA unique pointer to the DOF table created from the boundary mesh.
time_decay_parameterA parameter that defines the scaling factor of the nodal force. It should be a monotonic decrease parameter, meaning it should decrease over time.

Definition at line 21 of file ReleaseNodalForce.cpp.

26 : variable_id_(variable_id),
27 boundary_mesh_(boundary_mesh),
28 dof_table_(std::move(dof_table)),
29 time_decay_parameter_(time_decay_parameter)
30{
31}
ParameterLib::Parameter< double > const & time_decay_parameter_
MeshLib::Mesh const & boundary_mesh_
std::unique_ptr< NumLib::LocalToGlobalIndexMap > const dof_table_

Member Function Documentation

◆ applyNaturalBC()

void ProcessLib::ReleaseNodalForce::applyNaturalBC ( const double t,
std::vector< GlobalVector * > const & x,
int const process_id,
GlobalMatrix * K,
GlobalVector & b,
GlobalMatrix * Jac )
overridevirtual

Applies the released nodal force boundary condition. This method scales the nodal forces by the release parameter and applies them to the global RHS vector b.

Parameters
tThe current time.
xThe global solution vector (not used in this case).
process_idThe ID of the process.
KThe global stiffness matrix (not used in this case).
bThe global RHS vector to which the nodal forces are applied.
JacThe global Jacobian matrix (not used in this case).

Reimplemented from ProcessLib::BoundaryCondition.

Definition at line 80 of file ReleaseNodalForce.cpp.

85{
86 DBUG("Apply ReleaseNodalForce.");
87
89 {
90 return;
91 }
92
93 std::vector<double> release_values;
94 release_values.reserve(global_indices_.size());
95
96 for (std::size_t i = 0; i < global_indices_.size(); ++i)
97 {
99 auto const* node = boundary_nodes_[i];
100 pos.setNodeID(node->getID());
101 pos.setCoordinates(*node);
102
103 release_values.push_back(
104 -(1.0 - time_decay_parameter_(t, pos).front()) *
106 }
107
108 b.add(global_indices_, release_values);
109}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:77
void setNodeID(std::size_t node_id)
void setCoordinates(MathLib::Point3d const &coordinates)
std::vector< double > initial_release_nodal_force_
std::vector< GlobalIndexType > global_indices_
std::vector< MeshLib::Node const * > boundary_nodes_

References MathLib::EigenVector::add(), boundary_nodes_, DBUG(), global_indices_, initial_release_nodal_force_, ParameterLib::SpatialPosition::setCoordinates(), ParameterLib::SpatialPosition::setNodeID(), and time_decay_parameter_.

◆ set()

void ProcessLib::ReleaseNodalForce::set ( GlobalVector const * r_neq)

Definition at line 33 of file ReleaseNodalForce.cpp.

34{
35 auto const& number_of_components =
36 dof_table_->getNumberOfVariableComponents(variable_id_);
37
39 {
40 return;
41 }
42
44 number_of_components);
45
46 // Iterate over all nodes in the boundary mesh and set the initial
47 // released nodal.
48 for (auto const* node : boundary_mesh_.getNodes())
49 {
50 auto const node_id = node->getID();
53 for (int component_id = 0; component_id < number_of_components;
54 ++component_id)
55 {
56 auto const global_index =
57 dof_table_->getGlobalIndex(l, variable_id_, component_id);
58 if (global_index == NumLib::MeshComponentMap::nop)
59 {
60 continue;
61 }
62 // For the DDC approach (e.g. with PETSc option), the negative
63 // index of global_index means that the entry by that index is
64 // a ghost one, which should be dropped. Especially for PETSc
65 // routines MatZeroRows and MatZeroRowsColumns, which are
66 // called to apply the Dirichlet BC, the negative index is not
67 // accepted like other matrix or vector PETSc routines.
68 // Therefore, the following if-condition is applied.
69 if (global_index >= 0) [[likely]]
70 {
71 global_indices_.push_back(global_index);
73 r_neq->get(global_index));
74 boundary_nodes_.push_back(node);
75 }
76 }
77 }
78}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:108
std::size_t getID() const
Get id of the mesh.
Definition Mesh.h:123
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:102
static constexpr NUMLIB_EXPORT GlobalIndexType const nop

References boundary_mesh_, boundary_nodes_, dof_table_, MathLib::EigenVector::get(), MeshLib::Mesh::getID(), MeshLib::Mesh::getNodes(), MeshLib::Mesh::getNumberOfNodes(), global_indices_, initial_release_nodal_force_, MeshLib::Node, NumLib::MeshComponentMap::nop, and variable_id_.

Referenced by ProcessLib::BoundaryConditionCollection::setReleaseNodalForces().

Member Data Documentation

◆ boundary_mesh_

MeshLib::Mesh const& ProcessLib::ReleaseNodalForce::boundary_mesh_
private

Definition at line 143 of file ReleaseNodalForce.h.

Referenced by set().

◆ boundary_nodes_

std::vector<MeshLib::Node const*> ProcessLib::ReleaseNodalForce::boundary_nodes_
private

Definition at line 163 of file ReleaseNodalForce.h.

Referenced by applyNaturalBC(), and set().

◆ dof_table_

std::unique_ptr<NumLib::LocalToGlobalIndexMap> const ProcessLib::ReleaseNodalForce::dof_table_
private

Definition at line 145 of file ReleaseNodalForce.h.

Referenced by set().

◆ global_indices_

std::vector<GlobalIndexType> ProcessLib::ReleaseNodalForce::global_indices_
private

Definition at line 162 of file ReleaseNodalForce.h.

Referenced by applyNaturalBC(), and set().

◆ initial_release_nodal_force_

std::vector<double> ProcessLib::ReleaseNodalForce::initial_release_nodal_force_
private

Initial nodal force values at the boundary nodes. This vector is used to store the initial nodal forces before they are scaled by the release parameter.

Definition at line 160 of file ReleaseNodalForce.h.

Referenced by applyNaturalBC(), and set().

◆ time_decay_parameter_

ParameterLib::Parameter<double> const& ProcessLib::ReleaseNodalForce::time_decay_parameter_
private

A monotonically decreasing parameter that defines the scaling factor of the nodal force. It is a time- and position-dependent release parameter \( g(t, \mathbf{x}) \) representing the progressive removal of support over time, e.g., \( g(0, \mathbf{x}) = 1 \) and \( g(t_e, \mathbf{x}) = 0 \), where \( t_e \) is the end time of excavation, and \( \frac{\partial g}{\partial t} < 0 \).

Definition at line 155 of file ReleaseNodalForce.h.

Referenced by applyNaturalBC().

◆ variable_id_

int const ProcessLib::ReleaseNodalForce::variable_id_
private

Definition at line 141 of file ReleaseNodalForce.h.

Referenced by set().


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