26namespace ThermoMechanics
28template <
int DisplacementDim>
31 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
32 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
33 unsigned const integration_order,
34 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
38 bool const use_monolithic_scheme)
39 :
Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
40 integration_order, std::move(process_variables),
41 std::move(secondary_variables), use_monolithic_scheme),
42 _process_data(std::move(process_data))
51 std::make_unique<MeshLib::IntegrationPointWriter>(
58 std::make_unique<MeshLib::IntegrationPointWriter>(
65 std::make_unique<MeshLib::IntegrationPointWriter>(
72template <
int DisplacementDim>
78template <
int DisplacementDim>
81 const int process_id)
const
85 if (_use_monolithic_scheme ||
86 process_id == _process_data.mechanics_process_id)
88 auto const& l = *_local_to_global_index_map;
89 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
90 &l.getGhostIndices(), &this->_sparsity_pattern};
94 auto const& l = *_local_to_global_index_map_single_component;
95 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
96 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
101template <
int DisplacementDim>
107 if (_use_monolithic_scheme)
109 constructMonolithicProcessDofTable();
112 constructDofTableOfSpecifiedProcessStaggeredScheme(
113 _process_data.mechanics_process_id);
117 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
118 *_mesh_subset_all_nodes};
119 _local_to_global_index_map_single_component.reset(
121 std::move(all_mesh_subsets_single_component),
125 if (!_use_monolithic_scheme)
127 _sparsity_pattern_with_single_component =
129 *_local_to_global_index_map_single_component, _mesh);
133template <
int DisplacementDim>
137 unsigned const integration_order)
145 auto add_secondary_variable = [&](std::string
const& name,
146 int const num_components,
147 auto get_ip_values_function)
149 _secondary_variables.addSecondaryVariable(
153 std::move(get_ip_values_function)));
156 add_secondary_variable(
"sigma",
158 DisplacementDim>::RowsAtCompileTime,
159 &LocalAssemblerInterface::getIntPtSigma);
161 add_secondary_variable(
"epsilon",
163 DisplacementDim>::RowsAtCompileTime,
164 &LocalAssemblerInterface::getIntPtEpsilon);
171 add_secondary_variable);
179 *_local_to_global_index_map);
182template <
int DisplacementDim>
184 std::map<
int, std::shared_ptr<MaterialPropertyLib::Medium>>
const& media)
186 if (_use_monolithic_scheme)
188 const int process_id_of_thermomechanics = 0;
189 initializeProcessBoundaryConditionsAndSourceTerms(
190 *_local_to_global_index_map, process_id_of_thermomechanics, media);
196 initializeProcessBoundaryConditionsAndSourceTerms(
197 *_local_to_global_index_map_single_component,
198 _process_data.heat_conduction_process_id, media);
201 initializeProcessBoundaryConditionsAndSourceTerms(
202 *_local_to_global_index_map, _process_data.mechanics_process_id, media);
205template <
int DisplacementDim>
207 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
208 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
211 DBUG(
"Assemble ThermoMechanicsProcess.");
213 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
214 _local_to_global_index_map.get()};
219 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
223template <
int DisplacementDim>
226 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
227 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
230 DBUG(
"AssembleJacobian ThermoMechanicsProcess.");
232 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
234 if (_use_monolithic_scheme)
237 "Assemble the Jacobian of ThermoMechanics for the monolithic"
239 dof_tables.emplace_back(_local_to_global_index_map.get());
244 if (process_id == _process_data.heat_conduction_process_id)
247 "Assemble the Jacobian equations of heat conduction process in "
248 "ThermoMechanics for the staggered scheme.");
253 "Assemble the Jacobian equations of mechanical process in "
254 "ThermoMechanics for the staggered scheme.");
258 if (_process_data.heat_conduction_process_id ==
261 dof_tables.emplace_back(
262 _local_to_global_index_map_single_component.get());
263 dof_tables.emplace_back(_local_to_global_index_map.get());
267 dof_tables.emplace_back(_local_to_global_index_map.get());
268 dof_tables.emplace_back(
269 _local_to_global_index_map_single_component.get());
275 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
276 process_id, &b, &Jac);
279 auto copyRhs = [&](
int const variable_id,
auto& output_vector)
281 if (_use_monolithic_scheme)
283 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
285 std::negate<double>());
289 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
291 std::negate<double>());
294 if (_use_monolithic_scheme ||
295 process_id == _process_data.heat_conduction_process_id)
297 copyRhs(0, *_heat_flux);
299 if (_use_monolithic_scheme ||
300 process_id == _process_data.mechanics_process_id)
302 copyRhs(1, *_nodal_forces);
306template <
int DisplacementDim>
308 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
309 const int process_id)
311 DBUG(
"PreTimestep ThermoMechanicsProcess.");
313 assert(process_id < 2);
315 if (process_id == _process_data.mechanics_process_id)
319 getActiveElementIDs(), *_local_to_global_index_map, *x[process_id],
325template <
int DisplacementDim>
327 std::vector<GlobalVector*>
const& x,
328 std::vector<GlobalVector*>
const& x_prev,
double const t,
double const dt,
329 int const process_id)
336 DBUG(
"PostTimestep ThermoMechanicsProcess.");
337 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
338 auto const n_processes = x.size();
339 dof_tables.reserve(n_processes);
340 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
342 dof_tables.push_back(&getDOFTable(process_id));
347 getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
350template <
int DisplacementDim>
354 if (_process_data.mechanics_process_id == process_id)
356 return *_local_to_global_index_map;
360 return *_local_to_global_index_map_single_component;
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Global vector based on Eigen vector.
bool isAxiallySymmetric() const
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Properties & getProperties()
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
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)
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Handles configuration of several secondary variables from the project file.
Local assembler of ThermoMechanics process.
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
MeshLib::PropertyVector< double > * _nodal_forces
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
bool isLinear() const override
void constructDofTable() override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
ThermoMechanicsProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const ¶meters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, ThermoMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
MeshLib::PropertyVector< double > * _heat_flux
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *M, GlobalMatrix *K, GlobalVector *b)
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector *b, GlobalMatrix *Jac)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > getSigma() const =0
virtual std::vector< double > getEpsilonMechanical() const =0
virtual std::vector< double > getEpsilon() const =0