ExaDG
Loading...
Searching...
No Matches
material_handler.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_STRUCTURE_MATERIAL_MATERIAL_HANDLER_H_
23#define INCLUDE_EXADG_STRUCTURE_MATERIAL_MATERIAL_HANDLER_H_
24
25// deal.II
26#include <deal.II/matrix_free/matrix_free.h>
27
28// ExaDG
29#include <exadg/structure/material/library/st_venant_kirchhoff.h>
30#include <exadg/structure/material/material.h>
31#include <exadg/structure/user_interface/material_descriptor.h>
32
33namespace ExaDG
34{
35namespace Structure
36{
37template<int dim, typename Number>
39{
40public:
41 typedef std::pair<dealii::types::material_id, std::shared_ptr<Material<dim, Number>>> Pair;
42 typedef std::map<dealii::types::material_id, std::shared_ptr<Material<dim, Number>>> Materials;
43
44 MaterialHandler() : dof_index(0)
45 {
46 }
47
48 void
49 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
50 unsigned int const dof_index,
51 unsigned int const quad_index,
52 std::shared_ptr<MaterialDescriptor const> material_descriptor,
53 bool const large_deformation)
54 {
55 this->dof_index = dof_index;
56 this->material_descriptor = material_descriptor;
57
58 for(auto iter = material_descriptor->begin(); iter != material_descriptor->end(); ++iter)
59 {
60 dealii::types::material_id id = iter->first;
61 std::shared_ptr<MaterialData> data = iter->second;
62 MaterialType type = data->type;
63
64 switch(type)
65 {
66 case MaterialType::Undefined:
67 {
68 AssertThrow(false, dealii::ExcMessage("Material type is undefined."));
69 break;
70 }
71 case MaterialType::StVenantKirchhoff:
72 {
73 std::shared_ptr<StVenantKirchhoffData<dim>> data_svk =
74 std::static_pointer_cast<StVenantKirchhoffData<dim>>(data);
75 material_map.insert(
76 Pair(id,
78 matrix_free, dof_index, quad_index, *data_svk, large_deformation)));
79 break;
80 }
81 default:
82 {
83 AssertThrow(false, dealii::ExcMessage("Specified material type is not implemented."));
84 break;
85 }
86 }
87 }
88 }
89
90 void
91 reinit(dealii::MatrixFree<dim, Number> const & matrix_free, unsigned int const cell)
92 {
93 auto mid = matrix_free.get_cell_iterator(cell, 0, dof_index)->material_id();
94
95#ifdef DEBUG
96 for(unsigned int v = 1; v < matrix_free.n_active_entries_per_cell_batch(cell); v++)
97 AssertThrow(mid == matrix_free.get_cell_iterator(cell, v)->material_id(),
98 dealii::ExcMessage("You have to categorize cells according to their materials!"));
99#endif
100
101 material = material_map[mid];
102 }
103
104 std::shared_ptr<Material<dim, Number>>
105 get_material() const
106 {
107 return material;
108 }
109
110private:
111 unsigned int dof_index;
112
113 std::shared_ptr<MaterialDescriptor const> material_descriptor;
114 Materials material_map;
115
116 // pointer to material of current cell
117 std::shared_ptr<Material<dim, Number>> material;
118};
119
120} // namespace Structure
121} // namespace ExaDG
122
123#endif /* INCLUDE_EXADG_STRUCTURE_MATERIAL_MATERIAL_HANDLER_H_ */
Definition material_handler.h:39
Definition st_venant_kirchhoff.h:60
Definition driver.cpp:33