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