![]() |
OGS
|
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:
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.
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>
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 . | |
![]() | |
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< GlobalIndexType > | global_indices_ |
std::vector< MeshLib::Node const * > | boundary_nodes_ |
|
explicit |
Constructs a released nodal force boundary condition.
variable_id | The ID of the variable to which the nodal forces are applied. |
boundary_mesh | The mesh representing the boundary where the nodal forces are applied. |
dof_table | A unique pointer to the DOF table created from the boundary mesh. |
time_decay_parameter | A 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.
|
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
.
t | The current time. |
x | The global solution vector (not used in this case). |
process_id | The ID of the process. |
K | The global stiffness matrix (not used in this case). |
b | The global RHS vector to which the nodal forces are applied. |
Jac | The global Jacobian matrix (not used in this case). |
Reimplemented from ProcessLib::BoundaryCondition.
Definition at line 80 of file ReleaseNodalForce.cpp.
References MathLib::EigenVector::add(), boundary_nodes_, DBUG(), global_indices_, initial_release_nodal_force_, ParameterLib::SpatialPosition::setCoordinates(), ParameterLib::SpatialPosition::setNodeID(), and time_decay_parameter_.
void ProcessLib::ReleaseNodalForce::set | ( | GlobalVector const * | r_neq | ) |
Definition at line 33 of file ReleaseNodalForce.cpp.
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().
|
private |
Definition at line 143 of file ReleaseNodalForce.h.
Referenced by set().
|
private |
Definition at line 163 of file ReleaseNodalForce.h.
Referenced by applyNaturalBC(), and set().
|
private |
Definition at line 145 of file ReleaseNodalForce.h.
Referenced by set().
|
private |
Definition at line 162 of file ReleaseNodalForce.h.
Referenced by applyNaturalBC(), and set().
|
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().
|
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().
|
private |
Definition at line 141 of file ReleaseNodalForce.h.
Referenced by set().