SPEED
MAKE_STRAIN_TENSOR.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
63
64
65 subroutine make_strain_tensor(nn,ct,ww,dd,&
66 alfa11,alfa12,alfa13,alfa21,alfa22,alfa23,&
67 alfa31,alfa32,alfa33,beta11,beta12,beta13,&
68 beta21,beta22,beta23,beta31,beta32,beta33,&
69 gamma1,gamma2,gamma3,delta1,delta2,delta3,&
70 ux,uy,uz,&
71 duxdx,duydx,duzdx,&
72 duxdy,duydy,duzdy,&
73 duxdz,duydz,duzdz)
74
75 implicit none
76
77 integer*4 :: nn
78 integer*4 :: i,j,k,p,q,r
79
80 real*8 :: alfa11,alfa12,alfa13,alfa21,alfa22,alfa23,alfa31,alfa32,alfa33
81 real*8 :: beta11,beta12,beta13,beta21,beta22,beta23,beta31,beta32,beta33
82 real*8 :: gamma1,gamma2,gamma3,delta1,delta2,delta3
83 real*8 :: dxdx,dxdy,dxdz,dydx,dydy,dydz,dzdx,dzdy,dzdz,det_j
84 real*8 :: t1ux,t2ux,t3ux,t1uy,t2uy,t3uy,t1uz,t2uz,t3uz
85
86 real*8, dimension(nn) :: ct,ww
87 real*8, dimension(nn,nn) :: dd
88
89 real*8, dimension(nn,nn,nn) :: ux,uy,uz
90 real*8, dimension(nn,nn,nn) :: duxdx,duxdy,duxdz,duydx,duydy,duydz,duzdx,duzdy,duzdz
91
92! STRAIN CALCULATION
93
94 do r = 1,nn
95 do q = 1,nn
96 do p = 1,nn
97 t1ux = 0.d0
98 t1uy = 0.d0
99 t1uz = 0.d0
100 t2ux = 0.d0
101 t2uy = 0.d0
102 t2uz = 0.d0
103 t3ux = 0.d0
104 t3uy = 0.d0
105 t3uz = 0.d0
106
107 dxdx = alfa11 + beta12*ct(r) + beta13*ct(q) &
108 + gamma1*ct(q)*ct(r)
109 dydx = alfa21 + beta22*ct(r) + beta23*ct(q) &
110 + gamma2*ct(q)*ct(r)
111 dzdx = alfa31 + beta32*ct(r) + beta33*ct(q) &
112 + gamma3*ct(q)*ct(r)
113
114 dxdy = alfa12 + beta11*ct(r) + beta13*ct(p) &
115 + gamma1*ct(r)*ct(p)
116 dydy = alfa22 + beta21*ct(r) + beta23*ct(p) &
117 + gamma2*ct(r)*ct(p)
118 dzdy = alfa32 + beta31*ct(r) + beta33*ct(p) &
119 + gamma3*ct(r)*ct(p)
120
121 dxdz = alfa13 + beta11*ct(q) + beta12*ct(p) &
122 + gamma1*ct(p)*ct(q)
123 dydz = alfa23 + beta21*ct(q) + beta22*ct(p) &
124 + gamma2*ct(p)*ct(q)
125 dzdz = alfa33 + beta31*ct(q) + beta32*ct(p) &
126 + gamma3*ct(p)*ct(q)
127
128 det_j = dxdz * (dydx*dzdy - dzdx*dydy) &
129 - dydz * (dxdx*dzdy - dzdx*dxdy) &
130 + dzdz * (dxdx*dydy - dydx*dxdy)
131
132
133 do i = 1,nn
134 t1ux = t1ux + ux(i,q,r) * dd(p,i)
135 t1uy = t1uy + uy(i,q,r) * dd(p,i)
136 t1uz = t1uz + uz(i,q,r) * dd(p,i)
137 enddo
138
139 do j = 1,nn
140 t2ux = t2ux + ux(p,j,r) * dd(q,j)
141 t2uy = t2uy + uy(p,j,r) * dd(q,j)
142 t2uz = t2uz + uz(p,j,r) * dd(q,j)
143 enddo
144
145 do k = 1,nn
146 t3ux = t3ux + ux(p,q,k) * dd(r,k)
147 t3uy = t3uy + uy(p,q,k) * dd(r,k)
148 t3uz = t3uz + uz(p,q,k) * dd(r,k)
149 enddo
150
151
152 duxdx(p,q,r) = 1.0d0 / det_j *(&
153 ((dydy*dzdz - dydz*dzdy) * t1ux) + &
154 ((dydz*dzdx - dydx*dzdz) * t2ux) + &
155 ((dydx*dzdy - dydy*dzdx) * t3ux))
156
157 duydx(p,q,r) = 1.0d0 / det_j *(&
158 ((dydy*dzdz - dydz*dzdy) * t1uy) + &
159 ((dydz*dzdx - dydx*dzdz) * t2uy) + &
160 ((dydx*dzdy - dydy*dzdx) * t3uy))
161
162 duzdx(p,q,r) = 1.0d0 / det_j *(&
163 ((dydy*dzdz - dydz*dzdy) * t1uz) + &
164 ((dydz*dzdx - dydx*dzdz) * t2uz) + &
165 ((dydx*dzdy - dydy*dzdx) * t3uz))
166
167 duxdy(p,q,r) = 1.0d0 / det_j *(&
168 ((dzdy*dxdz - dzdz*dxdy) * t1ux) + &
169 ((dzdz*dxdx - dzdx*dxdz) * t2ux) + &
170 ((dzdx*dxdy - dzdy*dxdx) * t3ux))
171
172 duydy(p,q,r) = 1.0d0 / det_j *(&
173 ((dzdy*dxdz - dzdz*dxdy) * t1uy) + &
174 ((dzdz*dxdx - dzdx*dxdz) * t2uy) + &
175 ((dzdx*dxdy - dzdy*dxdx) * t3uy))
176
177 duzdy(p,q,r) = 1.0d0 / det_j *(&
178 ((dzdy*dxdz - dzdz*dxdy) * t1uz) + &
179 ((dzdz*dxdx - dzdx*dxdz) * t2uz) + &
180 ((dzdx*dxdy - dzdy*dxdx) * t3uz))
181
182 duxdz(p,q,r) = 1.0d0 / det_j *(&
183 ((dxdy*dydz - dxdz*dydy) * t1ux) + &
184 ((dxdz*dydx - dxdx*dydz) * t2ux) + &
185 ((dxdx*dydy - dxdy*dydx) * t3ux))
186
187 duydz(p,q,r) = 1.0d0 / det_j *(&
188 ((dxdy*dydz - dxdz*dydy) * t1uy) + &
189 ((dxdz*dydx - dxdx*dydz) * t2uy) + &
190 ((dxdx*dydy - dxdy*dydx) * t3uy))
191
192 duzdz(p,q,r) = 1.0d0 / det_j *(&
193 ((dxdy*dydz - dxdz*dydy) * t1uz) + &
194 ((dxdz*dydx - dxdx*dydz) * t2uz) + &
195 ((dxdx*dydy - dxdy*dydx) * t3uz))
196
197 enddo
198 enddo
199 enddo
200
201 return
202
203 end subroutine make_strain_tensor
204
subroutine make_strain_tensor(nn, ct, ww, dd, alfa11, alfa12, alfa13, alfa21, alfa22, alfa23
Computes spatial derivatives of displacement.