OGS
CentralDifferencesJacobianAssembler.cpp
Go to the documentation of this file.
1
12
13#include "BaseLib/ConfigTree.h"
14#include "BaseLib/Error.h"
17
18namespace ProcessLib
19{
21 std::vector<double>&& absolute_epsilons)
22 : _absolute_epsilons(std::move(absolute_epsilons))
23{
24 if (_absolute_epsilons.empty())
25 {
26 OGS_FATAL("No values for the absolute epsilons have been given.");
27 }
28}
29
31 LocalAssemblerInterface& local_assembler, const double t, double const dt,
32 const std::vector<double>& local_x_data,
33 const std::vector<double>& local_x_prev_data,
34 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
35{
36 std::vector<double> local_M_data(local_Jac_data.size());
37 std::vector<double> local_K_data(local_Jac_data.size());
38
39 // TODO do not check in every call.
40 if (local_x_data.size() % _absolute_epsilons.size() != 0)
41 {
43 "The number of specified epsilons ({:d}) and the number of local "
44 "d.o.f.s ({:d}) do not match, i.e., the latter is not divisible by "
45 "the former.",
46 _absolute_epsilons.size(), local_x_data.size());
47 }
48
49 auto const num_r_c =
50 static_cast<Eigen::MatrixXd::Index>(local_x_data.size());
51
52 auto const local_x =
53 MathLib::toVector<Eigen::VectorXd>(local_x_data, num_r_c);
54 auto const local_x_prev =
55 MathLib::toVector<Eigen::VectorXd>(local_x_prev_data, num_r_c);
56 Eigen::VectorXd const local_xdot = (local_x - local_x_prev) / dt;
57
58 auto local_Jac =
59 MathLib::createZeroedMatrix(local_Jac_data, num_r_c, num_r_c);
60 _local_x_perturbed_data = local_x_data;
61
62 auto const num_dofs_per_component =
63 local_x_data.size() / _absolute_epsilons.size();
64
65 auto const vds = local_assembler.getVectorDeformationSegment();
66
67 // Residual res := M xdot + K x - b
68 // Computing Jac := dres/dx
69 // = M dxdot/dx + dM/dx xdot + K dx/dx + dK/dx x - db/dx
70 // with dxdot/dx = 1/dt and dx/dx = 1
71 // (Note: dM/dx and dK/dx actually have the second and
72 // third index transposed.)
73 // The loop computes the dM/dx, dK/dx and db/dx terms, the rest is computed
74 // afterwards.
75 for (Eigen::MatrixXd::Index i = 0; i < num_r_c; ++i)
76 {
77 if ((vds != std::nullopt) && i >= vds->start_index &&
78 (i < (vds->start_index + vds->size)))
79 {
80 // Avoid to compute the analytic block
81 continue;
82 }
83
84 // assume that local_x_data is ordered by component.
85 auto const component = i / num_dofs_per_component;
86 auto const eps = _absolute_epsilons[component];
87
89 local_assembler.assemble(t, dt, _local_x_perturbed_data,
90 local_x_prev_data, local_M_data, local_K_data,
91 local_b_data);
92
93 _local_x_perturbed_data[i] = local_x_data[i] - eps;
94 local_assembler.assemble(t, dt, _local_x_perturbed_data,
95 local_x_prev_data, _local_M_data,
97
98 _local_x_perturbed_data[i] = local_x_data[i];
99
100 if (!local_M_data.empty())
101 {
102 auto const local_M_p =
103 MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
104 auto const local_M_m =
105 MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
106 local_Jac.col(i).noalias() +=
107 // dM/dxi * x_dot
108 (local_M_p - local_M_m) * local_xdot / (2.0 * eps);
109 local_M_data.clear();
110 _local_M_data.clear();
111 }
112 if (!local_K_data.empty())
113 {
114 auto const local_K_p =
115 MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
116 auto const local_K_m =
117 MathLib::toMatrix(_local_K_data, num_r_c, num_r_c);
118 local_Jac.col(i).noalias() +=
119 // dK/dxi * x
120 (local_K_p - local_K_m) * local_x / (2.0 * eps);
121 local_K_data.clear();
122 _local_K_data.clear();
123 }
124 if (!local_b_data.empty())
125 {
126 auto const local_b_p =
127 MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
128 auto const local_b_m =
130 local_Jac.col(i).noalias() -=
131 // db/dxi
132 (local_b_p - local_b_m) / (2.0 * eps);
133 local_b_data.clear();
134 _local_b_data.clear();
135 }
136 }
137
138 // Assemble with unperturbed local x, i.e. compute M dxdot/dx + K dx/dx =
139 // M/dt + K
140 local_assembler.assemble(t, dt, local_x_data, local_x_prev_data,
141 local_M_data, local_K_data, local_b_data);
142
143 // Compute remaining terms of the Jacobian.
144 if (!local_M_data.empty())
145 {
146 auto local_M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
147 local_Jac.noalias() += local_M / dt;
148 }
149 if (!local_K_data.empty())
150 {
151 auto local_K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
152 local_Jac.noalias() += local_K;
153 }
154
155 // Move the M and K contributions to the residuum for evaluation of nodal
156 // forces, flow rates, and the like. Cleaning up the M's and K's storage so
157 // it is not accounted for twice.
158 auto b = [&]()
159 {
160 if (!local_b_data.empty())
161 {
162 return MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
163 }
165 num_r_c);
166 }();
167
168 if (!local_M_data.empty())
169 {
170 auto M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
171 b -= M * local_xdot;
172 local_M_data.clear();
173 }
174 if (!local_K_data.empty())
175 {
176 auto K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
177
178 // Note: The deformation segment of \c b is already computed as
179 // int{B^T sigma}dA, which is identical to K_uu * u. Therefore the
180 // corresponding K block is set to zero.
181 if (vds != std::nullopt)
182 {
183 K.block(vds->start_index, vds->start_index, vds->size, vds->size)
184 .setZero();
185 }
186
187 b -= K * local_x;
188 local_K_data.clear();
189 }
190}
191
192std::unique_ptr<CentralDifferencesJacobianAssembler>
194{
196 config.checkConfigParameter("type", "CentralDifferences");
197
198 // TODO make non-optional.
200 auto rel_eps = config.getConfigParameterOptional<std::vector<double>>(
201 "relative_epsilons");
203 auto comp_mag = config.getConfigParameterOptional<std::vector<double>>(
204 "component_magnitudes");
205
206 if (!!rel_eps != !!comp_mag)
207 {
208 OGS_FATAL(
209 "Either both or none of <relative_epsilons> and "
210 "<component_magnitudes> have to be specified.");
211 }
212
213 std::vector<double> abs_eps;
214
215 if (rel_eps)
216 {
217 if (rel_eps->size() != comp_mag->size())
218 {
219 OGS_FATAL(
220 "The numbers of components of <relative_epsilons> and "
221 "<component_magnitudes> do not match.");
222 }
223
224 abs_eps.resize(rel_eps->size());
225 for (std::size_t i = 0; i < rel_eps->size(); ++i)
226 {
227 abs_eps[i] = (*rel_eps)[i] * (*comp_mag)[i];
228 }
229 }
230 else
231 {
232 // By default 1e-8 is used as epsilon for all components.
233 // TODO: remove this default value.
234 abs_eps.emplace_back(1e-8);
235 }
236
237 return std::make_unique<CentralDifferencesJacobianAssembler>(
238 std::move(abs_eps));
239}
240
241std::unique_ptr<AbstractJacobianAssembler>
243{
244 return std::make_unique<CentralDifferencesJacobianAssembler>(*this);
245}
246
247} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
std::optional< T > getConfigParameterOptional(std::string const &param) const
void checkConfigParameter(std::string const &param, std::string_view const value) const
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 std::optional< VectorSegment > getVectorDeformationSegment() 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)
std::unique_ptr< CentralDifferencesJacobianAssembler > createCentralDifferencesJacobianAssembler(BaseLib::ConfigTree const &config)