OGS
CentralDifferencesJacobianAssembler.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
7#include "BaseLib/Error.h"
10
11namespace ProcessLib
12{
14 std::vector<double>&& absolute_epsilons)
15 : AbstractJacobianAssembler(std::move(absolute_epsilons))
16{
17}
18
20 LocalAssemblerInterface& local_assembler, const double t, double const dt,
21 const std::vector<double>& local_x_data,
22 const std::vector<double>& local_x_prev_data,
23 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
24{
25 std::vector<double> local_M_data(local_Jac_data.size());
26 std::vector<double> local_K_data(local_Jac_data.size());
27
28 auto const num_r_c =
29 static_cast<Eigen::MatrixXd::Index>(local_x_data.size());
30
31 auto const local_x =
32 MathLib::toVector<Eigen::VectorXd>(local_x_data, num_r_c);
33 auto const local_x_prev =
34 MathLib::toVector<Eigen::VectorXd>(local_x_prev_data, num_r_c);
35 Eigen::VectorXd const local_xdot = (local_x - local_x_prev) / dt;
36
37 auto local_Jac =
38 MathLib::createZeroedMatrix(local_Jac_data, num_r_c, num_r_c);
39 _local_x_perturbed_data = local_x_data;
40
41 assert(this->non_deformation_component_ids_.size() > 0);
42 auto const num_dofs_per_component =
43 local_x_data.size() / this->non_deformation_component_ids_.size();
44
45 auto const nved = local_assembler.getNumberOfVectorElementsForDeformation();
46
47 // Residual res := M xdot + K x - b
48 // Computing Jac := dres/dx
49 // = M dxdot/dx + dM/dx xdot + K dx/dx + dK/dx x - db/dx
50 // with dxdot/dx = 1/dt and dx/dx = 1
51 // (Note: dM/dx and dK/dx actually have the second and
52 // third index transposed.)
53 // The loop computes the dM/dx, dK/dx and db/dx terms, the rest is computed
54 // afterwards. The loop skips the entries corresponding to the deformation
55 // part of the solution vector if a vector segment size is given by nved.
56 // This is to avoid recomputing the analytic block of the deformation
57 // process.
58 auto const num_purterbated_colums = num_r_c - nved;
59 auto const perturbations = this->getVariableComponentEpsilonsView();
60 for (Eigen::MatrixXd::Index i = 0; i < num_purterbated_colums; ++i)
61 {
62 // assume that local_x_data is ordered by component.
63 auto const component = i / num_dofs_per_component;
64 auto const eps = perturbations[component];
65
67 local_assembler.assemble(t, dt, _local_x_perturbed_data,
68 local_x_prev_data, local_M_data, local_K_data,
69 local_b_data);
70
71 _local_x_perturbed_data[i] = local_x_data[i] - eps;
72 local_assembler.assemble(t, dt, _local_x_perturbed_data,
73 local_x_prev_data, _local_M_data,
75
76 _local_x_perturbed_data[i] = local_x_data[i];
77
78 if (!local_M_data.empty())
79 {
80 auto const local_M_p =
81 MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
82 auto const local_M_m =
83 MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
84 local_Jac.col(i).noalias() +=
85 // dM/dxi * x_dot
86 (local_M_p - local_M_m) * local_xdot / (2.0 * eps);
87 local_M_data.clear();
88 _local_M_data.clear();
89 }
90 if (!local_K_data.empty())
91 {
92 auto const local_K_p =
93 MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
94 auto const local_K_m =
95 MathLib::toMatrix(_local_K_data, num_r_c, num_r_c);
96 local_Jac.col(i).noalias() +=
97 // dK/dxi * x
98 (local_K_p - local_K_m) * local_x / (2.0 * eps);
99 local_K_data.clear();
100 _local_K_data.clear();
101 }
102 if (!local_b_data.empty())
103 {
104 auto const local_b_p =
105 MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
106 auto const local_b_m =
108 local_Jac.col(i).noalias() -=
109 // db/dxi
110 (local_b_p - local_b_m) / (2.0 * eps);
111 local_b_data.clear();
112 _local_b_data.clear();
113 }
114 }
115
116 // Assemble with unperturbed local x, i.e. compute M dxdot/dx + K dx/dx =
117 // M/dt + K
118 local_assembler.assemble(t, dt, local_x_data, local_x_prev_data,
119 local_M_data, local_K_data, local_b_data);
120
121 // Compute remaining terms of the Jacobian.
122 if (!local_M_data.empty())
123 {
124 auto local_M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
125 local_Jac.noalias() += local_M / dt;
126 }
127 if (!local_K_data.empty())
128 {
129 auto local_K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
130 local_Jac.noalias() += local_K;
131 }
132
133 // Move the M and K contributions to the residuum for evaluation of nodal
134 // forces, flow rates, and the like. Cleaning up the M's and K's storage so
135 // it is not accounted for twice.
136 auto b = [&]()
137 {
138 if (!local_b_data.empty())
139 {
140 return MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
141 }
143 num_r_c);
144 }();
145
146 if (!local_M_data.empty())
147 {
148 auto M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
149 b -= M * local_xdot;
150 local_M_data.clear();
151 }
152 if (!local_K_data.empty())
153 {
154 auto K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
155
156 // Note: The deformation segment of \c b is already computed as
157 // int{B^T sigma}dA, which is identical to K_uu * u. Therefore the
158 // corresponding K block is set to zero.
159 if (nved != 0)
160 {
161 auto const dm_start_index = num_purterbated_colums;
162 auto const dm_size = nved;
163 K.block(dm_start_index, dm_start_index, dm_size, dm_size).setZero();
164 }
165
166 b -= K * local_x;
167 local_K_data.clear();
168 }
169}
170
171std::unique_ptr<AbstractJacobianAssembler>
173{
174 return std::make_unique<CentralDifferencesJacobianAssembler>(*this);
175}
176
177} // namespace ProcessLib
AbstractJacobianAssembler(std::vector< double > const &&absolute_epsilons)
CentralDifferencesJacobianAssembler(std::vector< double > &&absolute_epsilons)
std::unique_ptr< AbstractJacobianAssembler > copy() const override
void assembleWithJacobian(LocalAssemblerInterface &local_assembler, double const t, double const dt, std::vector< double > const &local_x_data, std::vector< double > const &local_x_prev_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
virtual int getNumberOfVectorElementsForDeformation() const
virtual void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)