1 #ifndef __COMMON_SIGMA_NUT_KERNEL_H__
2 #define __COMMON_SIGMA_NUT_KERNEL_H__
52 T * __restrict__ delta,
54 T * __restrict__ mult,
59 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
60 const int str = blockDim.x * gridDim.x;
62 for (
int i = idx;
i < n;
i += str) {
63 T sigG11, sigG12, sigG13, sigG22, sigG23, sigG33;
64 T sigma1, sigma2, sigma3;
65 T Invariant1, Invariant2, Invariant3;
66 T alpha1, alpha2, alpha3;
71 pi_3 = 4.0/3.0*atan(1.0);
76 g21[
i] = g21[
i] * mult[
i];
79 g31[
i] = g31[
i] * mult[
i];
80 g32[
i] = g32[
i] * mult[
i];
96 if (abs(sigG11) < eps) {
99 if (abs(sigG12) < eps) {
102 if (abs(sigG13) < eps) {
105 if (abs(sigG22) < eps) {
108 if (abs(sigG23) < eps) {
111 if (abs(sigG33) < eps) {
115 if (abs(sigG12*sigG12 + \
116 sigG13*sigG13 + sigG23*sigG23) < eps) {
117 sigma1 = sqrt(
max(
max(
max(sigG11, sigG22), sigG33), 0.0));
118 sigma3 = sqrt(
max(min(min(sigG11, sigG22), sigG33), 0.0));
119 Invariant1 = sigG11 + sigG22 + sigG33;
120 sigma2 = sqrt(abs(Invariant1 - sigma1*sigma1 - sigma3*sigma3));
122 Invariant1 = sigG11 + sigG22 + sigG33;
123 Invariant2 = sigG11*sigG22 + sigG11*sigG33 + sigG22*sigG33 - \
124 (sigG12*sigG12 + sigG13*sigG13 + sigG23*sigG23);
125 Invariant3 = sigG11*sigG22*sigG33 + \
126 2.0*sigG12*sigG13*sigG23 - \
127 (sigG11*sigG23*sigG23 + sigG22*sigG13*sigG13 + \
128 sigG33*sigG12*sigG12);
130 Invariant1 =
max(Invariant1, 0.0);
131 Invariant2 =
max(Invariant2, 0.0);
132 Invariant3 =
max(Invariant3, 0.0);
134 alpha1 = Invariant1*Invariant1/9.0 - Invariant2/3.0;
136 alpha1 =
max(alpha1, 0.0);
138 alpha2 = Invariant1*Invariant1*Invariant1/27.0 - \
139 Invariant1*Invariant2/6.0 + Invariant3/2.0;
141 tmp1 = alpha2/sqrt(alpha1 * alpha1 * alpha1);
144 sigma1 = sqrt(
max(Invariant1/3.0 + sqrt(alpha1), 0.0));
146 sigma3 = sqrt(Invariant1/3.0 - 2.0*sqrt(alpha1));
148 }
else if (tmp1 >= 1.0) {
149 sigma1 = sqrt(
max(Invariant1/3.0 + 2.0*sqrt(alpha1), 0.0));
150 sigma2 = sqrt(Invariant1/3.0 - sqrt(alpha1));
154 alpha3 = acos(tmp1)/3.0;
156 if (abs(Invariant3) < eps) {
157 sigma1 = sqrt(
max(Invariant1/3.0 + \
158 2.0*sqrt(alpha1)*cos(alpha3), 0.0));
159 sigma2 = sqrt(abs(Invariant1 - sigma1*sigma1));
162 sigma1 = sqrt(
max(Invariant1/3.0 + \
163 2.0*sqrt(alpha1)*cos(alpha3), 0.0));
164 sigma2 = sqrt(Invariant1/3.0 - \
165 2.0*sqrt(alpha1)*cos(pi_3 + alpha3));
166 sigma3 = sqrt(abs(Invariant1 - \
167 sigma1*sigma1-sigma2*sigma2));
173 Dsigma = sigma3*(sigma1 - sigma2)*(sigma2 - sigma3)/(sigma1*sigma1);
178 Dsigma =
max(Dsigma, 0.0);
180 nut[
i] = (c*delta[
i])*(c*delta[
i]) * Dsigma;
__global__ void sigma_nut_compute(T *__restrict__ g11, T *__restrict__ g12, T *__restrict__ g13, T *__restrict__ g21, T *__restrict__ g22, T *__restrict__ g23, T *__restrict__ g31, T *__restrict__ g32, T *__restrict__ g33, T *__restrict__ delta, T *__restrict__ nut, T *__restrict__ mult, const T c, const T eps, const int n)
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g23
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g22
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g13
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g12
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g33
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g11