OGS
CentralDifferencesJacobianAssembler.cpp
Go to the documentation of this file.
1 
12 
13 #include "BaseLib/ConfigTree.h"
14 #include "BaseLib/Error.h"
17 
18 namespace 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_xdot_data, const double dxdot_dx,
34  const double dx_dx, std::vector<double>& local_M_data,
35  std::vector<double>& local_K_data, std::vector<double>& local_b_data,
36  std::vector<double>& local_Jac_data)
37 {
38  // TODO do not check in every call.
39  if (local_x_data.size() % _absolute_epsilons.size() != 0)
40  {
41  OGS_FATAL(
42  "The number of specified epsilons ({:d}) and the number of local "
43  "d.o.f.s ({:d}) do not match, i.e., the latter is not divisible by "
44  "the former.",
45  _absolute_epsilons.size(), local_x_data.size());
46  }
47 
48  auto const num_r_c =
49  static_cast<Eigen::MatrixXd::Index>(local_x_data.size());
50 
51  auto const local_x =
52  MathLib::toVector<Eigen::VectorXd>(local_x_data, num_r_c);
53  auto const local_xdot =
54  MathLib::toVector<Eigen::VectorXd>(local_xdot_data, num_r_c);
55 
56  auto local_Jac =
57  MathLib::createZeroedMatrix(local_Jac_data, num_r_c, num_r_c);
58  _local_x_perturbed_data = local_x_data;
59  _local_xdot_perturbed_data = local_xdot_data;
60 
61  auto const num_dofs_per_component =
62  local_x_data.size() / _absolute_epsilons.size();
63 
64  // Residual res := M xdot + K x - b
65  // Computing Jac := dres/dx
66  // = M dxdot/dx + dM/dx xdot + K dx/dx + dK/dx x - db/dx
67  // (Note: dM/dx and dK/dx actually have the second and
68  // third index transposed.)
69  // The loop computes the dM/dx, dK/dx and db/dx terms, the rest is computed
70  // afterwards.
71  for (Eigen::MatrixXd::Index i = 0; i < num_r_c; ++i)
72  {
73  // assume that local_x_data is ordered by component.
74  auto const component = i / num_dofs_per_component;
75  auto const eps = _absolute_epsilons[component];
76 
77  _local_x_perturbed_data[i] += eps;
78  _local_xdot_perturbed_data[i] = local_xdot_data[i] + eps / dt;
79  local_assembler.assemble(t, dt, _local_x_perturbed_data,
80  _local_xdot_perturbed_data, local_M_data,
81  local_K_data, local_b_data);
82 
83  _local_x_perturbed_data[i] = local_x_data[i] - eps;
84  _local_xdot_perturbed_data[i] = local_xdot_data[i] - eps / dt;
85  local_assembler.assemble(t, dt, _local_x_perturbed_data,
88 
89  _local_x_perturbed_data[i] = local_x_data[i];
90  _local_xdot_perturbed_data[i] = local_xdot_data[i];
91 
92  if (!local_M_data.empty())
93  {
94  auto const local_M_p =
95  MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
96  auto const local_M_m =
97  MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
98  local_Jac.col(i).noalias() +=
99  // dM/dxi * x_dot
100  (local_M_p - local_M_m) * local_xdot / (2.0 * eps);
101  local_M_data.clear();
102  _local_M_data.clear();
103  }
104  if (!local_K_data.empty())
105  {
106  auto const local_K_p =
107  MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
108  auto const local_K_m =
109  MathLib::toMatrix(_local_K_data, num_r_c, num_r_c);
110  local_Jac.col(i).noalias() +=
111  // dK/dxi * x
112  (local_K_p - local_K_m) * local_x / (2.0 * eps);
113  local_K_data.clear();
114  _local_K_data.clear();
115  }
116  if (!local_b_data.empty())
117  {
118  auto const local_b_p =
119  MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
120  auto const local_b_m =
121  MathLib::toVector<Eigen::VectorXd>(_local_b_data, num_r_c);
122  local_Jac.col(i).noalias() -=
123  // db/dxi
124  (local_b_p - local_b_m) / (2.0 * eps);
125  local_b_data.clear();
126  _local_b_data.clear();
127  }
128  }
129 
130  // Assemble with unperturbed local x.
131  local_assembler.assemble(t, dt, local_x_data, local_xdot_data, local_M_data,
132  local_K_data, local_b_data);
133 
134  // Compute remaining terms of the Jacobian.
135  if (dxdot_dx != 0.0 && !local_M_data.empty())
136  {
137  auto local_M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
138  local_Jac.noalias() += local_M * dxdot_dx;
139  }
140  if (dx_dx != 0.0 && !local_K_data.empty())
141  {
142  auto local_K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
143  local_Jac.noalias() += local_K * dx_dx;
144  }
145 }
146 
147 std::unique_ptr<CentralDifferencesJacobianAssembler>
149 {
151  config.checkConfigParameter("type", "CentralDifferences");
152 
153  // TODO make non-optional.
155  auto rel_eps = config.getConfigParameterOptional<std::vector<double>>(
156  "relative_epsilons");
158  auto comp_mag = config.getConfigParameterOptional<std::vector<double>>(
159  "component_magnitudes");
160 
161  if (!!rel_eps != !!comp_mag)
162  {
163  OGS_FATAL(
164  "Either both or none of <relative_epsilons> and "
165  "<component_magnitudes> have to be specified.");
166  }
167 
168  std::vector<double> abs_eps;
169 
170  if (rel_eps)
171  {
172  if (rel_eps->size() != comp_mag->size())
173  {
174  OGS_FATAL(
175  "The numbers of components of <relative_epsilons> and "
176  "<component_magnitudes> do not match.");
177  }
178 
179  abs_eps.resize(rel_eps->size());
180  for (std::size_t i = 0; i < rel_eps->size(); ++i)
181  {
182  abs_eps[i] = (*rel_eps)[i] * (*comp_mag)[i];
183  }
184  }
185  else
186  {
187  // By default 1e-8 is used as epsilon for all components.
188  // TODO: remove this default value.
189  abs_eps.emplace_back(1e-8);
190  }
191 
192  return std::make_unique<CentralDifferencesJacobianAssembler>(
193  std::move(abs_eps));
194 }
195 
196 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
void checkConfigParameter(std::string const &param, T const &value) const
std::optional< T > getConfigParameterOptional(std::string const &param) const
CentralDifferencesJacobianAssembler(std::vector< double > &&absolute_epsilons)
void assembleWithJacobian(LocalAssemblerInterface &local_assembler, double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
virtual void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:72
std::unique_ptr< CentralDifferencesJacobianAssembler > createCentralDifferencesJacobianAssembler(BaseLib::ConfigTree const &config)