-
Notifications
You must be signed in to change notification settings - Fork 9
/
Copy pathstress.cu
140 lines (127 loc) · 4.24 KB
/
stress.cu
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
/**
* @author Christoph Schaefer [email protected]
*
* @section LICENSE
* Copyright (c) 2019 Christoph Schaefer
*
* This file is part of miluphcuda.
*
* miluphcuda is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* miluphcuda is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with miluphcuda. If not, see <http://www.gnu.org/licenses/>.
*
*/
#include "stress.h"
#include "parameter.h"
#include "miluph.h"
#include "timeintegration.h"
#include "linalg.h"
#if FRAGMENTATION
// if 1, then damage reduces the principal stresses
// if 0, then p<0 -> (1-d) p and S -> (1-d) S
// disabled for the time being
# define DAMAGE_ACTS_ON_PRINCIPAL_STRESSES 0
#else
# define DAMAGE_ACTS_ON_PRINCIPAL_STRESSES 0
#endif
#if DAMAGE_ACTS_ON_PRINCIPAL_STRESSES && ( MOHR_COULOMB_PLASTICITY || DRUCKER_PRAGER_PLASTICITY || COLLINS_PLASTICITY || COLLINS_PLASTICITY_SIMPLE )
# error Do not combine DAMAGE_ACTS_ON_PRINCIPAL_STRESSES and pressure-dependent yield strength models.
#endif
#if SOLID
// here we set the stress tensor sigma from pressure and deviatoric stress S
// note, that S was already lowered in plasticity
__global__ void set_stress_tensor(void)
{
register int i, inc, matId;
int d, e;
int niters;
double sigma[DIM][DIM];
# if DAMAGE_ACTS_ON_PRINCIPAL_STRESSES
double sigmatmp[DIM][DIM];
double rotation_matrix[DIM][DIM];
double main_stresses[DIM];
# endif
# if FRAGMENTATION
register double damage;
# endif
register double ptmp;
inc = blockDim.x * gridDim.x;
for (i = threadIdx.x + blockIdx.x * blockDim.x; i < numParticles; i += inc) {
matId = p_rhs.materialId[i];
niters = 0;
# if FRAGMENTATION
damage = p.damage_total[i];
if (damage > 1.0) damage = 1.0;
if (damage < 0.0) damage = 0.0;
# endif
# if DAMAGE_ACTS_ON_PRINCIPAL_STRESSES
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
sigmatmp[d][e] = 0.0;
sigma[d][e] = p.S[stressIndex(i, d, e)];
if (d == e) {
sigma[d][e] += -p.p[i];
}
}
}
// calculate main stresses
niters = calculate_all_eigenvalues(sigma, main_stresses, rotation_matrix);
for (d = 0; d < DIM; d++) {
sigmatmp[d][d] = main_stresses[d];
if (sigmatmp[d][d] > 0) {
sigmatmp[d][d] *= (1.0 - damage);
}
}
// rotate back the lowered principal stresses
multiply_matrix(sigmatmp, rotation_matrix, sigma);
transpose_matrix(rotation_matrix);
multiply_matrix(rotation_matrix, sigma, sigmatmp);
// sigmatmp now holds the stress tensor for particle i with damaged reduced stresses
copy_matrix(sigmatmp, sigma);
# else
// assemble stress tensor
# if COLLINS_PLASTICITY
// influence of damage on p < 0 already set in plasticity.cu
ptmp = p.p[i];
# elif FRAGMENTATION
if (p.p[i] < 0.0) {
// reduction of neg. pressure following Grady-Kipp model
ptmp = (1.0 - damage) * p.p[i];
} else {
ptmp = p.p[i];
}
# else
ptmp = p.p[i];
# endif
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
# if FRAGMENTATION && DAMAGE_ACTS_ON_S
// reduction of S following Grady-Kipp model
sigma[d][e] = (1.0 - damage) * p.S[stressIndex(i, d, e)];
# else
sigma[d][e] = p.S[stressIndex(i, d, e)];
# endif
if (d == e) { // the trace
sigma[d][e] -= ptmp;
}
}
}
# endif
// remember sigma
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
p_rhs.sigma[stressIndex(i,d,e)] = sigma[d][e];
}
}
}
}
#endif // SOLID