OGS 6.1.0-1721-g6382411ad
CentralDifferencesJacobianAssembler.cpp
Go to the documentation of this file.
1 
11 #include "BaseLib/Error.h"
14 
15 namespace ProcessLib
16 {
18  std::vector<double>&& absolute_epsilons)
19  : _absolute_epsilons(std::move(absolute_epsilons))
20 {
21  if (_absolute_epsilons.empty())
22  OGS_FATAL("No values for the absolute epsilons have been given.");
23 }
24 
26  LocalAssemblerInterface& local_assembler, const double t,
27  const std::vector<double>& local_x_data,
28  const std::vector<double>& local_xdot_data, const double dxdot_dx,
29  const double dx_dx, std::vector<double>& local_M_data,
30  std::vector<double>& local_K_data, std::vector<double>& local_b_data,
31  std::vector<double>& local_Jac_data)
32 {
33  // TODO do not check in every call.
34  if (local_x_data.size() % _absolute_epsilons.size() != 0) {
35  OGS_FATAL(
36  "The number of specified epsilons (%u) and the number of local "
37  "d.o.f.s (%u) do not match, i.e., the latter is not divisable by "
38  "the former.",
39  _absolute_epsilons.size(), local_x_data.size());
40  }
41 
42  auto const num_r_c =
43  static_cast<Eigen::MatrixXd::Index>(local_x_data.size());
44 
45  auto const local_x =
46  MathLib::toVector<Eigen::VectorXd>(local_x_data, num_r_c);
47  auto const local_xdot =
48  MathLib::toVector<Eigen::VectorXd>(local_xdot_data, num_r_c);
49 
50  auto local_Jac = MathLib::createZeroedMatrix(local_Jac_data,
51  num_r_c, num_r_c);
52  _local_x_perturbed_data = local_x_data;
53 
54  auto const num_dofs_per_component =
55  local_x_data.size() / _absolute_epsilons.size();
56 
57  // Residual res := M xdot + K x - b
58  // Computing Jac := dres/dx
59  // = M dxdot/dx + dM/dx xdot + K dx/dx + dK/dx x - db/dx
60  // (Note: dM/dx and dK/dx actually have the second and
61  // third index transposed.)
62  // The loop computes the dM/dx, dK/dx and db/dx terms, the rest is computed
63  // afterwards.
64  for (Eigen::MatrixXd::Index i = 0; i < num_r_c; ++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 = _absolute_epsilons[component];
69 
70  _local_x_perturbed_data[i] += eps;
71  local_assembler.assemble(t, _local_x_perturbed_data, local_M_data,
72  local_K_data, local_b_data);
73 
74  _local_x_perturbed_data[i] = local_x_data[i] - eps;
77 
78  _local_x_perturbed_data[i] = local_x_data[i];
79 
80  if (!local_M_data.empty()) {
81  auto const local_M_p =
82  MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
83  auto const local_M_m =
84  MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
85  local_Jac.col(i).noalias() +=
86  // dM/dxi * x_dot
87  (local_M_p - local_M_m) * local_xdot / (2.0 * eps);
88  local_M_data.clear();
89  _local_M_data.clear();
90  }
91  if (!local_K_data.empty()) {
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  auto const local_b_p =
104  MathLib::toVector<Eigen::VectorXd>(local_b_data, num_r_c);
105  auto const local_b_m =
106  MathLib::toVector<Eigen::VectorXd>(_local_b_data, num_r_c);
107  local_Jac.col(i).noalias() -=
108  // db/dxi
109  (local_b_p - local_b_m) / (2.0 * eps);
110  local_b_data.clear();
111  _local_b_data.clear();
112  }
113  }
114 
115  // Assemble with unperturbed local x.
116  local_assembler.assemble(t, local_x_data, local_M_data, local_K_data,
117  local_b_data);
118 
119  // Compute remaining terms of the Jacobian.
120  if (dxdot_dx != 0.0 && !local_M_data.empty()) {
121  auto local_M = MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
122  local_Jac.noalias() += local_M * dxdot_dx;
123  }
124  if (dx_dx != 0.0 && !local_K_data.empty()) {
125  auto local_K = MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
126  local_Jac.noalias() += local_K * dx_dx;
127  }
128 }
129 
130 std::unique_ptr<CentralDifferencesJacobianAssembler>
132 {
134  config.checkConfigParameter("type", "CentralDifferences");
135 
136  // TODO make non-optional.
138  auto rel_eps = config.getConfigParameterOptional<std::vector<double>>(
139  "relative_epsilons");
141  auto comp_mag = config.getConfigParameterOptional<std::vector<double>>(
142  "component_magnitudes");
143 
144  if (!!rel_eps != !!comp_mag) {
145  OGS_FATAL(
146  "Either both or none of <relative_epsilons> and "
147  "<component_magnitudes> have to be specified.");
148  }
149 
150  std::vector<double> abs_eps;
151 
152  if (rel_eps) {
153  if (rel_eps->size() != comp_mag->size()) {
154  OGS_FATAL(
155  "The numbers of components of <relative_epsilons> and "
156  "<component_magnitudes> do not match.");
157  }
158 
159  abs_eps.resize(rel_eps->size());
160  for (std::size_t i=0; i<rel_eps->size(); ++i) {
161  abs_eps[i] = (*rel_eps)[i] * (*comp_mag)[i];
162  }
163  } else {
164  // By default 1e-8 is used as epsilon for all components.
165  // TODO: remove this default value.
166  abs_eps.emplace_back(1e-8);
167  }
168 
169  return std::make_unique<CentralDifferencesJacobianAssembler>(
170  std::move(abs_eps));
171 }
172 
173 } // ProcessLib
std::unique_ptr< CentralDifferencesJacobianAssembler > createCentralDifferencesJacobianAssembler(BaseLib::ConfigTree const &config)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:31
void checkConfigParameter(std::string const &param, T const &value) const
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:71
boost::optional< T > getConfigParameterOptional(std::string const &param) const
#define OGS_FATAL(fmt,...)
Definition: Error.h:71
virtual void assemble(double const t, std::vector< double > const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void assembleWithJacobian(LocalAssemblerInterface &local_assembler, double const t, 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
CentralDifferencesJacobianAssembler(std::vector< double > &&absolute_epsilons)