OGS
LocalAssemblerInterface.cpp
Go to the documentation of this file.
1 
12 
13 #include <cassert>
14 
18 
19 namespace ProcessLib
20 {
22  double const /*t*/,
23  double const /*dt*/,
24  std::vector<double> const& /*local_x*/,
25  std::vector<double> const& /*local_xdot*/,
26  std::vector<double>& /*local_M_data*/,
27  std::vector<double>& /*local_K_data*/,
28  std::vector<double>& /*local_b_data*/)
29 {
30  OGS_FATAL(
31  "The assemble() function is not implemented in the local assembler.");
32 }
33 
35  double const /*t*/, double const /*dt*/, Eigen::VectorXd const& /*local_x*/,
36  Eigen::VectorXd const& /*local_xdot*/, int const /*process_id*/,
37  std::vector<double>& /*local_M_data*/,
38  std::vector<double>& /*local_K_data*/,
39  std::vector<double>& /*local_b_data*/)
40 {
41  OGS_FATAL(
42  "The assembleForStaggeredScheme() function is not implemented in the "
43  "local assembler.");
44 }
45 
47  double const /*t*/, double const /*dt*/,
48  std::vector<double> const& /*local_x*/,
49  std::vector<double> const& /*local_xdot*/, const double /*dxdot_dx*/,
50  const double /*dx_dx*/, std::vector<double>& /*local_M_data*/,
51  std::vector<double>& /*local_K_data*/,
52  std::vector<double>& /*local_b_data*/,
53  std::vector<double>& /*local_Jac_data*/)
54 {
55  OGS_FATAL(
56  "The assembleWithJacobian() function is not implemented in the local "
57  "assembler.");
58 }
59 
61  double const /*t*/, double const /*dt*/, Eigen::VectorXd const& /*local_x*/,
62  Eigen::VectorXd const& /*local_xdot*/, const double /*dxdot_dx*/,
63  const double /*dx_dx*/, int const /*process_id*/,
64  std::vector<double>& /*local_M_data*/,
65  std::vector<double>& /*local_K_data*/,
66  std::vector<double>& /*local_b_data*/,
67  std::vector<double>& /*local_Jac_data*/)
68 {
69  OGS_FATAL(
70  "The assembleWithJacobianForStaggeredScheme() function is not "
71  "implemented in the local assembler.");
72 }
73 
75  std::size_t const mesh_item_id,
76  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
77  double const t, double const dt, std::vector<GlobalVector*> const& x,
78  GlobalVector const& x_dot, int const process_id)
79 {
80  std::vector<double> local_x_vec;
81 
82  auto const n_processes = x.size();
83  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
84  {
85  auto const indices =
86  NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
87  assert(!indices.empty());
88  auto const local_solution = x[process_id]->get(indices);
89  local_x_vec.insert(std::end(local_x_vec), std::begin(local_solution),
90  std::end(local_solution));
91  }
92  auto const local_x = MathLib::toVector(local_x_vec);
93 
94  // Todo: A more decent way is to directly pass x_dots as done for x
95  auto const indices =
96  NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
97  auto const local_x_dot_vec = x_dot.get(indices);
98  auto const local_x_dot = MathLib::toVector(local_x_dot_vec);
99 
100  computeSecondaryVariableConcrete(t, dt, local_x, local_x_dot);
101 }
102 
104  std::size_t const mesh_item_id,
105  NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
106  double const t, bool const use_monolithic_scheme, int const process_id)
107 {
108  auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
109  auto const local_x = x.get(indices);
110 
111  setInitialConditionsConcrete(local_x, t, use_monolithic_scheme, process_id);
112 }
113 
115  std::size_t const /*mesh_item_id*/,
116  NumLib::LocalToGlobalIndexMap const& /*dof_table*/)
117 {
119 }
120 
122  std::size_t const mesh_item_id,
123  NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
124  double const t, double const delta_t)
125 {
126  auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
127  auto const local_x = x.get(indices);
128 
129  preTimestepConcrete(local_x, t, delta_t);
130 }
131 
133  std::size_t const mesh_item_id,
134  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
135  std::vector<GlobalVector*> const& x, double const t, double const dt)
136 {
137  std::vector<double> local_x_vec;
138 
139  auto const n_processes = x.size();
140  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
141  {
142  auto const indices =
143  NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
144  assert(!indices.empty());
145  auto const local_solution = x[process_id]->get(indices);
146  local_x_vec.insert(std::end(local_x_vec), std::begin(local_solution),
147  std::end(local_solution));
148  }
149  auto const local_x = MathLib::toVector(local_x_vec);
150 
151  postTimestepConcrete(local_x, t, dt);
152 }
153 
155  std::size_t const mesh_item_id,
156  NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
157  GlobalVector const& xdot, double const t, double const dt,
158  bool const use_monolithic_scheme, int const process_id)
159 {
160  auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
161  auto const local_x = x.get(indices);
162  auto const local_xdot = xdot.get(indices);
163 
164  postNonLinearSolverConcrete(local_x, local_xdot, t, dt,
165  use_monolithic_scheme, process_id);
166 }
167 
168 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
Global vector based on Eigen vector.
Definition: EigenVector.h:26
double get(IndexType rowId) const
get entry
Definition: EigenVector.h:62
virtual void preTimestepConcrete(std::vector< double > const &, double const, double const)
virtual void assembleForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void setInitialConditionsConcrete(std::vector< double > const &, double const, bool const, int const)
virtual void assembleWithJacobian(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)
virtual void computeSecondaryVariableConcrete(double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &)
virtual void postNonLinearSolverConcrete(std::vector< double > const &, std::vector< double > const &, double const, double const, bool const, int const)
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
void setInitialConditions(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id)
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)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
void postNonLinearSolver(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
virtual void postTimestepConcrete(Eigen::VectorXd const &, double const, double const)
virtual void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)