OGS
ReleaseNodalForce.cpp
Go to the documentation of this file.
1
12#include "ReleaseNodalForce.h"
13
15#include "MeshLib/Mesh.h"
18
19namespace ProcessLib
20{
22 int const variable_id,
23 MeshLib::Mesh const& boundary_mesh,
24 std::unique_ptr<NumLib::LocalToGlobalIndexMap>& dof_table,
25 ParameterLib::Parameter<double> const& time_decay_parameter)
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}
32
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}
79
81 std::vector<GlobalVector*> const& /*x*/,
82 int const /*process_id*/,
83 GlobalMatrix* /*K*/, GlobalVector& b,
84 GlobalMatrix* /*Jac*/)
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}
110
111} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of the Mesh class.
Global vector based on Eigen vector.
Definition EigenVector.h:26
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:77
double get(IndexType rowId) const
get entry
Definition EigenVector.h:59
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
void setNodeID(std::size_t node_id)
void setCoordinates(MathLib::Point3d const &coordinates)
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 relea...
std::vector< double > initial_release_nodal_force_
ParameterLib::Parameter< double > const & time_decay_parameter_
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.
std::vector< GlobalIndexType > global_indices_
MeshLib::Mesh const & boundary_mesh_
std::vector< MeshLib::Node const * > boundary_nodes_
void set(GlobalVector const *r_neq)
std::unique_ptr< NumLib::LocalToGlobalIndexMap > const dof_table_