forked from christophmschaefer/miluphcuda
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathartificial_stress.cu
109 lines (92 loc) · 3.15 KB
/
artificial_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
/**
* @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 "artificial_stress.h"
#include "parameter.h"
#include "miluph.h"
#include "linalg.h"
// transform to principal axes before applying artificial stress
#define PRINCIPAL_AXES_ARTIFICIAL_STRESS 1
#if ARTIFICIAL_STRESS
__global__ void compute_artificial_stress(int *interactions)
{
int i, inc;
int d, e;
int niters = 0;
int matId;
double max_ev = -1e300;
// the diagonalized tensors
double main_stresses[DIM];
double rotation_matrix[DIM][DIM];
double R[DIM][DIM];
double sigma[DIM][DIM];
double Rtmp[DIM][DIM];
inc = blockDim.x * gridDim.x;
for (i = threadIdx.x + blockIdx.x * blockDim.x; i < numParticles; i += inc) {
matId = p_rhs.materialId[i];
// build stress tensor from deviatoric stress and pressure
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
sigma[d][e] = p_rhs.sigma[stressIndex(i,d,e)];
}
}
#if PRINCIPAL_AXES_ARTIFICIAL_STRESS
niters = calculate_all_eigenvalues(sigma, main_stresses, rotation_matrix);
// determine the maximum stress
max_ev = main_stresses[0];
for (e = 1; e < DIM; e++) {
if (main_stresses[e] > max_ev) {
max_ev = main_stresses[e];
}
}
// now calculate the artificial stress from the main stresses
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
R[d][e] = 0.0;
}
if (main_stresses[d] > 0) {
R[d][d] = -matepsilon_stress[matId]*main_stresses[d];
}
}
// convert back in the original coordinate system with the rotation matrix
multiply_matrix(R, rotation_matrix, Rtmp);
transpose_matrix(rotation_matrix);
multiply_matrix(rotation_matrix, Rtmp, R);
// now save R for the particle
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
p_rhs.R[stressIndex(i,d,e)]= R[d][e];
}
}
#else
// no transformation, just reduce the standard stress
for (d = 0; d < DIM; d++) {
for (e = 0; e < DIM; e++) {
if (sigma[d][e] > 0) {
p_rhs.R[stressIndex(i,d,e)] = -matepsilon_stress[matId] * sigma[d][e];
}
}
}
#endif // PRINCIPAL_AXES_ARTIFICIAL_STRESS
}
}
#endif