OGS
SmallDeformationLocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <utility>
14#include <vector>
15
16#include "BaseLib/Error.h"
19
20namespace ProcessLib
21{
22namespace LIE
23{
24namespace SmallDeformation
25{
29{
30public:
33 std::size_t n_local_size, std::vector<unsigned> dofIndex_to_localIndex)
34 : _dofIndex_to_localIndex(std::move(dofIndex_to_localIndex))
35 {
36 _local_u.resize(n_local_size);
37 _local_b.resize(_local_u.size());
38 _local_J.resize(_local_u.size(), _local_u.size());
39 }
40
41 void assembleWithJacobian(double const t, double const dt,
42 std::vector<double> const& local_x_,
43 std::vector<double> const& /*local_x_prev*/,
44 std::vector<double>& /*local_M_data*/,
45 std::vector<double>& /*local_K_data*/,
46 std::vector<double>& local_b_data,
47 std::vector<double>& local_Jac_data) override
48 {
49 auto const local_dof_size = local_x_.size();
50
51 _local_u.setZero();
52 for (unsigned i = 0; i < local_dof_size; i++)
53 {
54 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
55 }
56 _local_b.setZero();
57 _local_J.setZero();
58
60
61 local_b_data.resize(local_dof_size);
62 for (unsigned i = 0; i < local_dof_size; i++)
63 {
64 local_b_data[i] = _local_b[_dofIndex_to_localIndex[i]];
65 }
66
67 local_Jac_data.resize(local_dof_size * local_dof_size);
68 for (unsigned i = 0; i < local_dof_size; i++)
69 {
70 for (unsigned j = 0; j < local_dof_size; j++)
71 {
72 local_Jac_data[i * local_dof_size + j] = _local_J(
74 }
75 }
76 }
77
78 virtual void assembleWithJacobian(double const /*t*/, double const /*dt*/,
79 Eigen::VectorXd const& /*local_u*/,
80 Eigen::VectorXd& /*local_b*/,
81 Eigen::MatrixXd& /*local_J*/)
82 {
84 "SmallDeformationLocalAssemblerInterface::assembleWithJacobian() "
85 "is not implemented");
86 }
87
89 double const t, double const /*dt*/, Eigen::VectorXd const& local_x,
90 Eigen::VectorXd const& /*local_x_prev*/) override
91 {
92 if (!_dofIndex_to_localIndex.empty())
93 {
94 _local_u.setZero();
95 for (auto i = 0; i < local_x.rows(); i++)
96 {
97 _local_u[_dofIndex_to_localIndex[i]] = local_x[i];
98 }
99 }
100
102 }
103
104 virtual std::vector<double> const& getIntPtSigma(
105 const double t,
106 std::vector<GlobalVector*> const& x,
107 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
108 std::vector<double>& cache) const = 0;
109
110 virtual std::vector<double> const& getIntPtEpsilon(
111 const double t,
112 std::vector<GlobalVector*> const& x,
113 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
114 std::vector<double>& cache) const = 0;
115
116 virtual std::vector<double> const& getIntPtFractureStress(
117 const double t,
118 std::vector<GlobalVector*> const& x,
119 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
120 std::vector<double>& cache) const = 0;
121
122 virtual std::vector<double> const& getIntPtFractureAperture(
123 const double t,
124 std::vector<GlobalVector*> const& x,
125 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
126 std::vector<double>& cache) const = 0;
127
128protected:
130 double const t, Eigen::VectorXd const& local_u) = 0;
131
132private:
133 Eigen::VectorXd _local_u;
134 Eigen::VectorXd _local_b;
135 Eigen::MatrixXd _local_J;
136 std::vector<unsigned> const _dofIndex_to_localIndex;
137};
138
139} // namespace SmallDeformation
140} // namespace LIE
141} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtFractureAperture(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
SmallDeformationLocalAssemblerInterface(std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
virtual void computeSecondaryVariableConcreteWithVector(double const t, Eigen::VectorXd const &local_u)=0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
void computeSecondaryVariableConcrete(double const t, double const, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
virtual void assembleWithJacobian(double const, double const, Eigen::VectorXd const &, Eigen::VectorXd &, Eigen::MatrixXd &)
virtual std::vector< double > const & getIntPtFractureStress(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0