OGS
ForwardDifferencesJacobianAssembler.cpp
Go to the documentation of this file.
1
11
12#include "BaseLib/Error.h"
15
16namespace ProcessLib
17{
19 std::vector<double>&& absolute_epsilons)
20 : _absolute_epsilons(std::move(absolute_epsilons))
21{
22 if (_absolute_epsilons.empty())
23 {
24 OGS_FATAL("No values for the absolute epsilons have been given.");
25 }
26}
27
29 LocalAssemblerInterface& local_assembler, const double t, double const dt,
30 const std::vector<double>& local_x_data,
31 const std::vector<double>& local_x_prev_data,
32 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
33{
34 std::vector<double> local_M_data(local_Jac_data.size());
35 std::vector<double> local_K_data(local_Jac_data.size());
36
37 // TODO do not check in every call.
38 if (local_x_data.size() % _absolute_epsilons.size() != 0)
39 {
41 "The number of specified epsilons ({:d}) and the number of local "
42 "d.o.f.s ({:d}) do not match, i.e., the latter is not divisible by "
43 "the former.",
44 _absolute_epsilons.size(), local_x_data.size());
45 }
46
47 auto const num_r_c =
48 static_cast<Eigen::MatrixXd::Index>(local_x_data.size());
49
50 auto const x = MathLib::toVector<Eigen::VectorXd>(local_x_data, num_r_c);
51 auto const x_prev =
52 MathLib::toVector<Eigen::VectorXd>(local_x_prev_data, num_r_c);
53
54 auto local_Jac =
55 MathLib::createZeroedMatrix(local_Jac_data, num_r_c, num_r_c);
56
57 auto const num_dofs_per_component =
58 local_x_data.size() / _absolute_epsilons.size();
59
60 // Assemble with unperturbed local x to get M0, K0, and b0 used in the
61 // finite differences below.
62 local_assembler.assemble(t, dt, local_x_data, local_x_prev_data,
63 local_M_data, local_K_data, local_b_data);
64
65 // Residual res := M xdot + K x - b
66 // Computing Jac := dres/dx
67 // = d(M xdot)/dx + d(K x)/dx - db/dx
68 for (Eigen::MatrixXd::Index i = 0; i < num_r_c; ++i)
69 {
70 // assume that local_x_data is ordered by component.
71 auto const component = i / num_dofs_per_component;
72 auto const eps = _absolute_epsilons[component];
73
74 // Assemble with perturbed local x.
75 _local_x_perturbed_data = local_x_data;
76 _local_x_perturbed_data[i] = local_x_data[i] + eps;
79
80 local_assembler.assemble(t, dt, _local_x_perturbed_data,
81 local_x_prev_data, _local_M_data,
83
84 if (!local_M_data.empty() && !_local_M_data.empty())
85 {
86 auto const local_M_0 =
87 MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
88 auto const local_M_p =
89 MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
90 local_Jac.col(i).noalias() +=
91 (local_M_p * (x_p - x_prev) - local_M_0 * (x - x_prev)) /
92 (dt * eps);
93 _local_M_data.clear();
94 }
95 if (!local_K_data.empty() && !_local_K_data.empty())
96 {
97 auto const local_K_0 =
98 MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
99 auto const local_K_p =
100 MathLib::toMatrix(_local_K_data, num_r_c, num_r_c);
101 local_Jac.col(i).noalias() +=
102 (local_K_p * x_p - local_K_0 * x) / eps;
103 _local_K_data.clear();
104 }
105 if (!local_b_data.empty() && !_local_b_data.empty())
106 {
107 auto const local_b_0 =
108 MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
109 auto const local_b_p =
111 local_Jac.col(i).noalias() -= (local_b_p - local_b_0) / eps;
112 _local_b_data.clear();
113 }
114 }
115
116 local_M_data.clear();
117 local_K_data.clear();
118 local_b_data.clear();
119 // Assemble with unperturbed local x to keep the internal assembler state
120 // for next iteration or time step.
121 local_assembler.assemble(t, dt, local_x_data, local_x_prev_data,
122 local_M_data, local_K_data, local_b_data);
123
124 // Move the M and K contributions to the residuum for evaluation of nodal
125 // forces, flow rates, and the like. Cleaning up the M's and K's storage so
126 // it is not accounted for twice.
127 auto b = [&]()
128 {
129 if (!local_b_data.empty())
130 {
131 return MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
132 }
134 num_r_c);
135 }();
136
137 if (!local_M_data.empty())
138 {
139 auto M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
140 b -= M * (x - x_prev) / dt;
141 local_M_data.clear();
142 }
143 if (!local_K_data.empty())
144 {
145 auto K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
146 b -= K * x;
147 local_K_data.clear();
148 }
149}
150
151std::unique_ptr<AbstractJacobianAssembler>
153{
154 return std::make_unique<ForwardDifferencesJacobianAssembler>(*this);
155}
156
157} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
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 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)