13 #include <Eigen/Dense>
20 class LocalToGlobalIndexMap;
25 struct CoupledSolutionsForStaggeredScheme;
26 struct LocalCoupledSolutions;
41 bool const use_monolithic_scheme,
42 int const process_id);
44 virtual void initialize(std::size_t
const mesh_item_id,
48 std::vector<double>
const& ){};
50 virtual void assemble(
double const t,
double const dt,
51 std::vector<double>
const& local_x,
52 std::vector<double>
const& local_xdot,
53 std::vector<double>& local_M_data,
54 std::vector<double>& local_K_data,
55 std::vector<double>& local_b_data);
58 Eigen::VectorXd
const& local_x,
59 Eigen::VectorXd
const& local_xdot,
61 std::vector<double>& local_M_data,
62 std::vector<double>& local_K_data,
63 std::vector<double>& local_b_data);
66 std::vector<double>
const& local_x,
67 std::vector<double>
const& local_xdot,
68 const double dxdot_dx,
const double dx_dx,
69 std::vector<double>& local_M_data,
70 std::vector<double>& local_K_data,
71 std::vector<double>& local_b_data,
72 std::vector<double>& local_Jac_data);
75 double const t,
double const dt, Eigen::VectorXd
const& local_x,
76 Eigen::VectorXd
const& local_xdot,
const double dxdot_dx,
77 const double dx_dx,
int const process_id,
78 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
79 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
82 std::size_t
const mesh_item_id,
83 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
84 double const t,
double const dt, std::vector<GlobalVector*>
const& x,
87 virtual void preTimestep(std::size_t
const mesh_item_id,
90 double const delta_t);
93 std::size_t
const mesh_item_id,
94 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
95 std::vector<GlobalVector*>
const& x,
double const t,
double const dt);
100 double const t,
double const dt,
101 bool const use_monolithic_scheme,
102 int const process_id);
105 std::vector<double>
const& )
116 std::vector<double>
const& )
const
118 return Eigen::Vector3d{};
125 std::vector<std::vector<double>>
const& )
const
127 return Eigen::Vector3d{};
132 std::vector<double>
const& ,
double const ,
133 bool const ,
int const )
140 double const ,
double const )
145 double const ,
double const )
150 std::vector<double>
const& ,
151 std::vector<double>
const& ,
double const ,
152 double const ,
bool const ,
160 Eigen::VectorXd
const& ,
161 Eigen::VectorXd
const& )
166 double const , std::vector<std::vector<double>>
const&
Definition of the Point3d class.
Global vector based on Eigen vector.
virtual void preTimestepConcrete(std::vector< double > const &, double const, double const)
virtual ~LocalAssemblerInterface()=default
virtual Eigen::Vector3d getFlux(MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
Fits to staggered scheme.
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 initializeConcrete()
virtual void computeSecondaryVariableConcrete(double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &)
virtual Eigen::Vector3d getFlux(MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual void preAssemble(double const, double const, std::vector< double > 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)
virtual void computeSecondaryVariableWithCoupledProcessConcrete(double const, std::vector< std::vector< double >> const &)
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints(std::vector< double > const &)