Plugins
Loading...
Searching...
No Matches
voigtian_notation.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-2020: G-CSC, Goethe University Frankfurt
3 * Author: Lukas Larisch
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef VOIGTIAN_H_
34#define VOIGTIAN_H_
35
36namespace ug{
37namespace SmallStrainMechanics{
38
39template <typename TDomain>
41public:
42 static const int dim = IMaterialLaw<TDomain>::dim;
43
45
46 void set_orthotropic(const number C11, const number C12, const number C13,
47 const number C22, const number C23, const number C33,
48 const number C44, const number C55, const number C66)
49 {
50 mat(0, 0) = C11; mat(0, 1) = C12; mat(0, 2) = C13;
51 mat(1, 0) = C12; mat(1, 1) = C22; mat(1, 2) = C23;
52 mat(2, 0) = C13; mat(2, 1) = C23; mat(2, 2) = C33;
53
54 c44 = C44;
55 c55 = C55;
56 c66 = C66;
57 }
58
59 //counterclock-wise
62
63 rot(0, 0) = 0; //cos 90
64 rot(0, 1) = -1; //-sin 90
65 rot(1, 0) = 1; //sin 90
66 rot(1, 1) = 0; //cos 90
67 rot(2, 2) = 1; //fixed axis
68
69 mat = rot * mat;
70
71 c44 = -c44;
72 number tmp = c55;
73 c55 = c66;
74 c66 = tmp;
75 }
76
78 {
79 // setze alle wWerte auf 0
80 for (size_t i = 0; i < (size_t) dim; ++i)
81 for (size_t j = 0; j < (size_t) dim; ++j)
82 for (size_t k = 0; k < (size_t) dim; ++k)
83 for (size_t l = 0; l < (size_t) dim; ++l)
84 elastTensorFunct[i][j][k][l] = 0.0;
85
86 elastTensorFunct[0][0][0][0] = mat(0, 0);
87 elastTensorFunct[0][0][1][1] = mat(0, 1);
88 elastTensorFunct[0][0][2][2] = mat(0, 2);
89 //elastTensorFunct[0][0][1][2] = mat(0, 3);
90 //elastTensorFunct[0][0][2][0] = mat(0, 4);
91 //elastTensorFunct[0][0][0][1] = mat(0, 5);
92
93 elastTensorFunct[1][1][0][0] = mat(1, 0);
94 elastTensorFunct[1][1][1][1] = mat(1, 1);
95 elastTensorFunct[1][1][2][2] = mat(1, 2);
96 //elastTensorFunct[1][1][1][2] = mat(1, 3);
97 //elastTensorFunct[1][1][2][0] = mat(1, 4);
98 //elastTensorFunct[1][1][0][1] = mat(1, 5);
99
100 elastTensorFunct[2][2][0][0] = mat(2, 0);
101 elastTensorFunct[2][2][1][1] = mat(2, 1);
102 elastTensorFunct[2][2][2][2] = mat(2, 2);
103 //elastTensorFunct[2][2][1][2] = mat(2, 3);
104 //elastTensorFunct[2][2][2][0] = mat(2, 4);
105 //elastTensorFunct[2][2][0][1] = mat(2, 5);
106
107 //elastTensorFunct[1][2][0][0] = mat(3, 0);
108 //elastTensorFunct[1][2][1][1] = mat(3, 1);
109 //elastTensorFunct[1][2][2][2] = mat(3, 2);
111 elastTensorFunct[1][2][1][2] = c44;
112 //elastTensorFunct[1][2][2][0] = mat(3, 4);
113 //elastTensorFunct[1][2][0][1] = mat(3, 5);
114
115 //elastTensorFunct[2][0][0][0] = mat(4, 0);
116 //elastTensorFunct[2][0][1][1] = mat(4, 1);
117 //elastTensorFunct[2][0][2][2] = mat(4, 2);
118 //elastTensorFunct[2][0][1][2] = mat(4, 3);
120 elastTensorFunct[2][0][2][0] = c55;
121 //elastTensorFunct[2][0][0][1] = mat(4, 5);
122
123 //elastTensorFunct[0][1][0][0] = mat(5, 0);
124 //elastTensorFunct[0][1][1][1] = mat(5, 1);
125 //elastTensorFunct[0][1][2][2] = mat(5, 2);
126 //elastTensorFunct[0][1][1][2] = mat(5, 3);
127 //elastTensorFunct[0][1][2][0] = mat(5, 4);
129 elastTensorFunct[0][1][0][1] = c66;
130
131 //symmetries: c_ijkl = c_jikl = c_ijlk
132 for(size_t i = 0; i < dim; ++i){
133 for(size_t j = 0; j < dim; ++j){
134 for(size_t k = 0; k < dim; ++k){
135 for(size_t l = 0; l < dim; ++l){
136 if(elastTensorFunct[i][j][k][l] != 0){
137 elastTensorFunct[j][i][k][l] = elastTensorFunct[i][j][k][l];
138 elastTensorFunct[j][i][l][k] = elastTensorFunct[i][j][k][l];
139 elastTensorFunct[i][j][l][k] = elastTensorFunct[i][j][k][l];
140 }
141 }
142 }
143 }
144 }
145 }
146
147private:
149 number c44; //1/G_xy
150 number c55; //1/G_yz
151 number c66; //1/G_xz
152};
153
154
155}// end of namespace SmallStrainMechanics
156}// end of namespace ug
157
158#endif /* Voigtian_H_ */
Definition mat_law_interface.h:44
Definition voigtian_notation.h:40
number c55
Definition voigtian_notation.h:150
number c44
Definition voigtian_notation.h:149
number c66
Definition voigtian_notation.h:151
void rotate_z_90_deg()
Definition voigtian_notation.h:60
void set_orthotropic(const number C11, const number C12, const number C13, const number C22, const number C23, const number C33, const number C44, const number C55, const number C66)
Definition voigtian_notation.h:46
DenseMatrix< FixedArray2< double, dim, dim > > mat
Definition voigtian_notation.h:148
void copy_to_tensor(MathTensor4< dim, dim, dim, dim > &elastTensorFunct)
Definition voigtian_notation.h:77
VoigtianMatrix()
Definition voigtian_notation.h:44
static const int dim
Definition voigtian_notation.h:42
double number
static void rot(number A[3][3], const number s, const number tau, const int i, const int j, const int k, const int l)