OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
s8zfint_reg.F File Reference
#include "implicit_f.inc"
#include "parit_c.inc"
#include "scr02_c.inc"
#include "scr18_c.inc"
#include "vectorize.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine s8zfint_reg (nloc_dmg, var_reg, nel, offg, vol, nc1, nc2, nc3, nc4, nc5, nc6, nc7, nc8, px1, px2, px3, px4, px5, px6, px7, px8, py1, py2, py3, py4, py5, py6, py7, py8, pz1, pz2, pz3, pz4, pz5, pz6, pz7, pz8, imat, h, wi, ip, itask, dt2t, vol0, nft)

Function/Subroutine Documentation

◆ s8zfint_reg()

subroutine s8zfint_reg ( type(nlocal_str_), target nloc_dmg,
intent(in) var_reg,
integer nel,
intent(in) offg,
intent(in) vol,
integer, dimension(nel) nc1,
integer, dimension(nel) nc2,
integer, dimension(nel) nc3,
integer, dimension(nel) nc4,
integer, dimension(nel) nc5,
integer, dimension(nel) nc6,
integer, dimension(nel) nc7,
integer, dimension(nel) nc8,
intent(in) px1,
intent(in) px2,
intent(in) px3,
intent(in) px4,
intent(in) px5,
intent(in) px6,
intent(in) px7,
intent(in) px8,
intent(in) py1,
intent(in) py2,
intent(in) py3,
intent(in) py4,
intent(in) py5,
intent(in) py6,
intent(in) py7,
intent(in) py8,
intent(in) pz1,
intent(in) pz2,
intent(in) pz3,
intent(in) pz4,
intent(in) pz5,
intent(in) pz6,
intent(in) pz7,
intent(in) pz8,
integer imat,
dimension(8), intent(in) h,
wi,
integer ip,
integer itask,
intent(inout) dt2t,
intent(in) vol0,
integer, intent(in) nft )

Definition at line 32 of file s8zfint_reg.F.

45C-----------------------------------------------
46C M o d u l e s
47C-----------------------------------------------
49C-----------------------------------------------
50C I m p l i c i t T y p e s
51C-----------------------------------------------
52#include "implicit_f.inc"
53C-----------------------------------------------
54C G l o b a l P a r a m e t e r s
55C-----------------------------------------------
56#include "parit_c.inc"
57#include "scr02_c.inc"
58#include "scr18_c.inc"
59C-----------------------------------------------
60C D u m m y A r g u m e n t s
61C-----------------------------------------------
62 INTEGER, INTENT(IN) :: NFT
63 INTEGER :: NEL,IMAT,IP,ITASK
64 INTEGER, DIMENSION(NEL) :: NC1,NC2,NC3,NC4,NC5,NC6,NC7,NC8
65 my_real, DIMENSION(NEL), INTENT(IN) ::
66 . offg
67 my_real, INTENT(INOUT) ::
68 . dt2t
69 my_real, DIMENSION(NEL), INTENT(IN) ::
70 . var_reg,px1,px2,px3,px4,px5,px6,px7,px8,
71 . py1,py2,py3,py4,py5,py6,py7,py8,pz1,pz2,
72 . pz3,pz4,pz5,pz6,pz7,pz8,vol,h(8),vol0
73 TYPE(NLOCAL_STR_), TARGET :: NLOC_DMG
75 . wi
76C-----------------------------------------------
77C L o c a l V a r i a b l e s
78C-----------------------------------------------
79 INTEGER I,II,K,NNOD,N1,N2,N3,N4,N5,N6,N7,N8,L_NLOC
80 my_real
81 . l2,xi,
82 . b1,b2,b3,b4,b5,b6,b7,b8,
83 . a1,a2,a3,a4,a5,a6,a7,a8,c1,c2,c3,c4,c5,c6,c7,c8,
84 . zeta,sspnl,dtnl,le_max,maxstif,minmasscal
85 my_real, DIMENSION(NEL) ::
86 . f1,f2,f3,f4,f5,f6,f7,f8,lc
87 my_real, DIMENSION(:) ,ALLOCATABLE ::
88 . btb11,btb12,btb13,btb14,btb15,btb16,btb17,btb18,
89 . btb22,btb23,btb24,btb25,btb26,btb27,btb28,btb33,
90 . btb34,btb35,btb36,btb37,btb38,btb44,btb45,btb46,
91 . btb47,btb48,btb55,btb56,btb57,btb58,btb66,btb67,
92 . btb68,btb77,btb78,btb88,sti1,sti2,sti3,sti4,sti5,
93 . sti6,sti7,sti8
94 INTEGER, DIMENSION(:), ALLOCATABLE ::
95 . POS1,POS2,POS3,POS4,POS5,POS6,POS7,POS8
96 my_real, POINTER, DIMENSION(:) ::
97 . vnl,fnl,unl,stifnl,mass,mass0,vnl0
98C=======================================================================
99 NULLIFY(mass)
100 NULLIFY(mass0)
101 l2 = nloc_dmg%LEN(imat)**2
102 xi = nloc_dmg%DAMP(imat)
103 nnod = nloc_dmg%NNOD
104 l_nloc = nloc_dmg%L_NLOC
105 zeta = nloc_dmg%DENS(imat)
106 sspnl = nloc_dmg%SSPNL(imat)
107 le_max = nloc_dmg%LE_MAX(imat) ! Maximal length of convergence
108 lc(1:nel) = zero
109 vnl => nloc_dmg%VNL(1:l_nloc)
110 vnl0 => nloc_dmg%VNL_OLD(1:l_nloc)
111 unl => nloc_dmg%UNL(1:l_nloc)
112 ALLOCATE(btb11(nel),btb12(nel),btb13(nel),btb14(nel),btb15(nel),
113 . btb16(nel),btb17(nel),btb18(nel),btb22(nel),btb23(nel),btb24(nel),
114 . btb25(nel),btb26(nel),btb27(nel),btb28(nel),btb33(nel),btb34(nel),
115 . btb35(nel),btb36(nel),btb37(nel),btb38(nel),btb44(nel),btb45(nel),
116 . btb46(nel),btb47(nel),btb48(nel),btb55(nel),btb56(nel),btb57(nel),
117 . btb58(nel),btb66(nel),btb67(nel),btb68(nel),btb77(nel),btb78(nel),
118 . btb88(nel),pos1(nel),pos2(nel),pos3(nel),pos4(nel),pos5(nel),
119 . pos6(nel),pos7(nel),pos8(nel))
120 ! For nodal timestep
121 IF (nodadt > 0) THEN
122 ! Non-local nodal stifness
123 ALLOCATE(sti1(nel),sti2(nel),sti3(nel),sti4(nel),sti5(nel),sti6(nel),
124 . sti7(nel),sti8(nel))
125 ! Non-local mass
126 mass => nloc_dmg%MASS(1:l_nloc)
127 ! Initial non-local mass
128 mass0 => nloc_dmg%MASS0(1:l_nloc)
129 ENDIF
130c
131 !-----------------------------------------------------------------------
132 ! Computation of the element volume and the BtB matrix product
133 !-----------------------------------------------------------------------
134 ! Loop over elements
135# include "vectorize.inc"
136 DO i=1,nel
137c
138 ! Number of the element nodes
139 n1 = nloc_dmg%IDXI(nc1(i))
140 n2 = nloc_dmg%IDXI(nc2(i))
141 n3 = nloc_dmg%IDXI(nc3(i))
142 n4 = nloc_dmg%IDXI(nc4(i))
143 n5 = nloc_dmg%IDXI(nc5(i))
144 n6 = nloc_dmg%IDXI(nc6(i))
145 n7 = nloc_dmg%IDXI(nc7(i))
146 n8 = nloc_dmg%IDXI(nc8(i))
147c
148 ! Recovering the position of the non-local d.o.f.s
149 pos1(i) = nloc_dmg%POSI(n1)
150 pos2(i) = nloc_dmg%POSI(n2)
151 pos3(i) = nloc_dmg%POSI(n3)
152 pos4(i) = nloc_dmg%POSI(n4)
153 pos5(i) = nloc_dmg%POSI(n5)
154 pos6(i) = nloc_dmg%POSI(n6)
155 pos7(i) = nloc_dmg%POSI(n7)
156 pos8(i) = nloc_dmg%POSI(n8)
157c
158 ! Computation of the product BtxB
159 btb11(i) = px1(i)**2 + py1(i)**2 + pz1(i)**2
160 btb12(i) = px1(i)*px2(i) + py1(i)*py2(i) + pz1(i)*pz2(i)
161 btb13(i) = px1(i)*px3(i) + py1(i)*py3(i) + pz1(i)*pz3(i)
162 btb14(i) = px1(i)*px4(i) + py1(i)*py4(i) + pz1(i)*pz4(i)
163 btb15(i) = px1(i)*px5(i) + py1(i)*py5(i) + pz1(i)*pz5(i)
164 btb16(i) = px1(i)*px6(i) + py1(i)*py6(i) + pz1(i)*pz6(i)
165 btb17(i) = px1(i)*px7(i) + py1(i)*py7(i) + pz1(i)*pz7(i)
166 btb18(i) = px1(i)*px8(i) + py1(i)*py8(i) + pz1(i)*pz8(i)
167 btb22(i) = px2(i)**2 + py2(i)**2 + pz2(i)**2
168 btb23(i) = px2(i)*px3(i) + py2(i)*py3(i) + pz2(i)*pz3(i)
169 btb24(i) = px2(i)*px4(i) + py2(i)*py4(i) + pz2(i)*pz4(i)
170 btb25(i) = px2(i)*px5(i) + py2(i)*py5(i) + pz2(i)*pz5(i)
171 btb26(i) = px2(i)*px6(i) + py2(i)*py6(i) + pz2(i)*pz6(i)
172 btb27(i) = px2(i)*px7(i) + py2(i)*py7(i) + pz2(i)*pz7(i)
173 btb28(i) = px2(i)*px8(i) + py2(i)*py8(i) + pz2(i)*pz8(i)
174 btb33(i) = px3(i)**2 + py3(i)**2 + pz3(i)**2
175 btb34(i) = px3(i)*px4(i) + py3(i)*py4(i) + pz3(i)*pz4(i)
176 btb35(i) = px3(i)*px5(i) + py3(i)*py5(i) + pz3(i)*pz5(i)
177 btb36(i) = px3(i)*px6(i) + py3(i)*py6(i) + pz3(i)*pz6(i)
178 btb37(i) = px3(i)*px7(i) + py3(i)*py7(i) + pz3(i)*pz7(i)
179 btb38(i) = px3(i)*px8(i) + py3(i)*py8(i) + pz3(i)*pz8(i)
180 btb44(i) = px4(i)**2 + py4(i)**2 + pz4(i)**2
181 btb45(i) = px4(i)*px5(i) + py4(i)*py5(i) + pz4(i)*pz5(i)
182 btb46(i) = px4(i)*px6(i) + py4(i)*py6(i) + pz4(i)*pz6(i)
183 btb47(i) = px4(i)*px7(i) + py4(i)*py7(i) + pz4(i)*pz7(i)
184 btb48(i) = px4(i)*px8(i) + py4(i)*py8(i) + pz4(i)*pz8(i)
185 btb55(i) = px5(i)**2 + py5(i)**2 + pz5(i)**2
186 btb56(i) = px5(i)*px6(i) + py5(i)*py6(i) + pz5(i)*pz6(i)
187 btb57(i) = px5(i)*px7(i) + py5(i)*py7(i) + pz5(i)*pz7(i)
188 btb58(i) = px5(i)*px8(i) + py5(i)*py8(i) + pz5(i)*pz8(i)
189 btb66(i) = px6(i)**2 + py6(i)**2 + pz6(i)**2
190 btb67(i) = px6(i)*px7(i) + py6(i)*py7(i) + pz6(i)*pz7(i)
191 btb68(i) = px6(i)*px8(i) + py6(i)*py8(i) + pz6(i)*pz8(i)
192 btb77(i) = px7(i)**2 + py7(i)**2 + pz7(i)**2
193 btb78(i) = px7(i)*px8(i) + py7(i)*py8(i) + pz7(i)*pz8(i)
194 btb88(i) = px8(i)**2 + py8(i)**2 + pz8(i)**2
195c
196 ENDDO
197c
198 !-----------------------------------------------------------------------
199 ! Computation of non-local forces
200 !-----------------------------------------------------------------------
201 ! Loop over elements
202# include "vectorize.inc"
203 DO i=1,nel
204c
205 ! If the element is not broken, normal computation
206 IF (offg(i)/=zero) THEN
207 ! Computation of LEN**2*BTB*Unl product and NTN*UNL, NTN*VNL
208 a1 = vol(i) * (h(1)*h(1)*unl(pos1(i)) + h(1)*h(2)*unl(pos2(i)) + h(1)*h(3)*unl(pos3(i))
209 . + h(1)*h(4)*unl(pos4(i)) + h(1)*h(5)*unl(pos5(i)) + h(1)*h(6)*unl(pos6(i))
210 . + h(1)*h(7)*unl(pos7(i)) + h(1)*h(8)*unl(pos8(i)))
211c
212 IF (nodadt == 0) THEN
213 a1 = a1 + vol(i) * xi * (h(1)*h(1)*vnl(pos1(i)) + h(1)*h(2)*vnl(pos2(i)) + h(1)*h(3)*vnl(pos3(i))
214 . + h(1)*h(4)*vnl(pos4(i)) + h(1)*h(5)*vnl(pos5(i)) + h(1)*h(6)*vnl(pos6(i))
215 . + h(1)*h(7)*vnl(pos7(i)) + h(1)*h(8)*vnl(pos8(i)))
216 ELSE
217 minmasscal = min(sqrt(mass(pos1(i))/mass0(pos1(i))),
218 . sqrt(mass(pos2(i))/mass0(pos2(i))),
219 . sqrt(mass(pos3(i))/mass0(pos3(i))),
220 . sqrt(mass(pos4(i))/mass0(pos4(i))),
221 . sqrt(mass(pos5(i))/mass0(pos5(i))),
222 . sqrt(mass(pos6(i))/mass0(pos6(i))),
223 . sqrt(mass(pos7(i))/mass0(pos7(i))),
224 . sqrt(mass(pos8(i))/mass0(pos8(i))))
225 a1 = a1 + vol(i) * xi * (h(1)*h(1)*minmasscal*vnl(pos1(i)) +
226 . h(1)*h(2)*minmasscal*vnl(pos2(i)) +
227 . h(1)*h(3)*minmasscal*vnl(pos3(i)) +
228 . h(1)*h(4)*minmasscal*vnl(pos4(i)) +
229 . h(1)*h(5)*minmasscal*vnl(pos5(i)) +
230 . h(1)*h(6)*minmasscal*vnl(pos6(i)) +
231 . h(1)*h(7)*minmasscal*vnl(pos7(i)) +
232 . h(1)*h(8)*minmasscal*vnl(pos8(i)))
233 ENDIF
234c
235 b1 = l2 * vol(i) * ( btb11(i)*unl(pos1(i)) + btb12(i)*unl(pos2(i))
236 . + btb13(i)*unl(pos3(i)) + btb14(i)*unl(pos4(i)) + btb15(i)*unl(pos5(i))
237 . + btb16(i)*unl(pos6(i)) + btb17(i)*unl(pos7(i)) + btb18(i)*unl(pos8(i)))
238c
239 c1 = vol(i) * h(1) * var_reg(i)
240c
241 a2 = vol(i) * (h(2)*h(1)*unl(pos1(i)) + h(2)*h(2)*unl(pos2(i)) + h(2)*h(3)*unl(pos3(i))
242 . + h(2)*h(4)*unl(pos4(i)) + h(2)*h(5)*unl(pos5(i)) + h(2)*h(6)*unl(pos6(i))
243 . + h(2)*h(7)*unl(pos7(i)) + h(2)*h(8)*unl(pos8(i)))
244c
245 IF (nodadt == 0) THEN
246 a2 = a2 + vol(i) * xi * (h(2)*h(1)*vnl(pos1(i)) + h(2)*h(2)*vnl(pos2(i)) + h(2)*h(3)*vnl(pos3(i))
247 . + h(2)*h(4)*vnl(pos4(i)) + h(2)*h(5)*vnl(pos5(i)) + h(2)*h(6)*vnl(pos6(i))
248 . + h(2)*h(7)*vnl(pos7(i)) + h(2)*h(8)*vnl(pos8(i)))
249 ELSE
250 a2 = a2 + vol(i) * xi * (h(2)*h(1)*minmasscal*vnl(pos1(i)) +
251 . h(2)*h(2)*minmasscal*vnl(pos2(i)) +
252 . h(2)*h(3)*minmasscal*vnl(pos3(i)) +
253 . h(2)*h(4)*minmasscal*vnl(pos4(i)) +
254 . h(2)*h(5)*minmasscal*vnl(pos5(i)) +
255 . h(2)*h(6)*minmasscal*vnl(pos6(i)) +
256 . h(2)*h(7)*minmasscal*vnl(pos7(i)) +
257 . h(2)*h(8)*minmasscal*vnl(pos8(i)))
258 ENDIF
259c
260 b2 = l2 * vol(i) * ( btb12(i)*unl(pos1(i)) + btb22(i)*unl(pos2(i))
261 . + btb23(i)*unl(pos3(i)) + btb24(i)*unl(pos4(i)) + btb25(i)*unl(pos5(i))
262 . + btb26(i)*unl(pos6(i)) + btb27(i)*unl(pos7(i)) + btb28(i)*unl(pos8(i)))
263c
264 c2 = vol(i) * h(2) * var_reg(i)
265c
266 a3 = vol(i) * (h(3)*h(1)*unl(pos1(i)) + h(3)*h(2)*unl(pos2(i)) + h(3)*h(3)*unl(pos3(i))
267 . + h(3)*h(4)*unl(pos4(i)) + h(3)*h(5)*unl(pos5(i)) + h(3)*h(6)*unl(pos6(i))
268 . + h(3)*h(7)*unl(pos7(i)) + h(3)*h(8)*unl(pos8(i)))
269c
270 IF (nodadt == 0) THEN
271 a3 = a3 + vol(i) * xi * (h(3)*h(1)*vnl(pos1(i)) + h(3)*h(2)*vnl(pos2(i)) + h(3)*h(3)*vnl(pos3(i))
272 . + h(3)*h(4)*vnl(pos4(i)) + h(3)*h(5)*vnl(pos5(i)) + h(3)*h(6)*vnl(pos6(i))
273 . + h(3)*h(7)*vnl(pos7(i)) + h(3)*h(8)*vnl(pos8(i)))
274 ELSE
275 a3 = a3 + vol(i) * xi * (h(3)*h(1)*minmasscal*vnl(pos1(i)) +
276 . h(3)*h(2)*minmasscal*vnl(pos2(i)) +
277 . h(3)*h(3)*minmasscal*vnl(pos3(i)) +
278 . h(3)*h(4)*minmasscal*vnl(pos4(i)) +
279 . h(3)*h(5)*minmasscal*vnl(pos5(i)) +
280 . h(3)*h(6)*minmasscal*vnl(pos6(i)) +
281 . h(3)*h(7)*minmasscal*vnl(pos7(i)) +
282 . h(3)*h(8)*minmasscal*vnl(pos8(i)))
283 ENDIF
284c
285 b3 = l2 * vol(i) * ( btb13(i)*unl(pos1(i)) + btb23(i)*unl(pos2(i))
286 . + btb33(i)*unl(pos3(i)) + btb34(i)*unl(pos4(i)) + btb35(i)*unl(pos5(i))
287 . + btb36(i)*unl(pos6(i)) + btb37(i)*unl(pos7(i)) + btb38(i)*unl(pos8(i)))
288c
289 c3 = vol(i) * h(3) * var_reg(i)
290c
291 a4 = vol(i) * (h(4)*h(1)*unl(pos1(i)) + h(4)*h(2)*unl(pos2(i)) + h(4)*h(3)*unl(pos3(i))
292 . + h(4)*h(4)*unl(pos4(i)) + h(4)*h(5)*unl(pos5(i)) + h(4)*h(6)*unl(pos6(i))
293 . + h(4)*h(7)*unl(pos7(i)) + h(4)*h(8)*unl(pos8(i)))
294c
295 IF (nodadt == 0) THEN
296 a4 = a4 + vol(i) * xi * (h(4)*h(1)*vnl(pos1(i)) + h(4)*h(2)*vnl(pos2(i)) + h(4)*h(3)*vnl(pos3(i))
297 . + h(4)*h(4)*vnl(pos4(i)) + h(4)*h(5)*vnl(pos5(i)) + h(4)*h(6)*vnl(pos6(i))
298 . + h(4)*h(7)*vnl(pos7(i)) + h(4)*h(8)*vnl(pos8(i)))
299 ELSE
300 a4 = a4 + vol(i) * xi * (h(4)*h(1)*minmasscal*vnl(pos1(i)) +
301 . h(4)*h(2)*minmasscal*vnl(pos2(i)) +
302 . h(4)*h(3)*minmasscal*vnl(pos3(i)) +
303 . h(4)*h(4)*minmasscal*vnl(pos4(i)) +
304 . h(4)*h(5)*minmasscal*vnl(pos5(i)) +
305 . h(4)*h(6)*minmasscal*vnl(pos6(i)) +
306 . h(4)*h(7)*minmasscal*vnl(pos7(i)) +
307 . h(4)*h(8)*minmasscal*vnl(pos8(i)))
308 ENDIF
309c
310 b4 = l2 * vol(i) * ( btb14(i)*unl(pos1(i)) + btb24(i)*unl(pos2(i))
311 . + btb34(i)*unl(pos3(i)) + btb44(i)*unl(pos4(i)) + btb45(i)*unl(pos5(i))
312 . + btb46(i)*unl(pos6(i)) + btb47(i)*unl(pos7(i)) + btb48(i)*unl(pos8(i)))
313c
314 c4 = vol(i) * h(4) * var_reg(i)
315c
316 a5 = vol(i) * (h(5)*h(1)*unl(pos1(i)) + h(5)*h(2)*unl(pos2(i)) + h(5)*h(3)*unl(pos3(i))
317 . + h(5)*h(4)*unl(pos4(i)) + h(5)*h(5)*unl(pos5(i)) + h(5)*h(6)*unl(pos6(i))
318 . + h(5)*h(7)*unl(pos7(i)) + h(5)*h(8)*unl(pos8(i)))
319c
320 IF (nodadt == 0) THEN
321 a5 = a5 + vol(i) * xi * (h(5)*h(1)*vnl(pos1(i)) + h(5)*h(2)*vnl(pos2(i)) + h(5)*h(3)*vnl(pos3(i))
322 . + h(5)*h(4)*vnl(pos4(i)) + h(5)*h(5)*vnl(pos5(i)) + h(5)*h(6)*vnl(pos6(i))
323 . + h(5)*h(7)*vnl(pos7(i)) + h(5)*h(8)*vnl(pos8(i)))
324 ELSE
325 a5 = a5 + vol(i) * xi * (h(5)*h(1)*minmasscal*vnl(pos1(i)) +
326 . h(5)*h(2)*minmasscal*vnl(pos2(i)) +
327 . h(5)*h(3)*minmasscal*vnl(pos3(i)) +
328 . h(5)*h(4)*minmasscal*vnl(pos4(i)) +
329 . h(5)*h(5)*minmasscal*vnl(pos5(i)) +
330 . h(5)*h(6)*minmasscal*vnl(pos6(i)) +
331 . h(5)*h(7)*minmasscal*vnl(pos7(i)) +
332 . h(5)*h(8)*minmasscal*vnl(pos8(i)))
333 ENDIF
334c
335 b5 = l2 * vol(i) * ( btb15(i)*unl(pos1(i)) + btb25(i)*unl(pos2(i))
336 . + btb35(i)*unl(pos3(i)) + btb45(i)*unl(pos4(i)) + btb55(i)*unl(pos5(i))
337 . + btb56(i)*unl(pos6(i)) + btb57(i)*unl(pos7(i)) + btb58(i)*unl(pos8(i)))
338c
339 c5 = vol(i) * h(5) * var_reg(i)
340c
341 a6 = vol(i) * (h(6)*h(1)*unl(pos1(i)) + h(6)*h(2)*unl(pos2(i)) + h(6)*h(3)*unl(pos3(i))
342 . + h(6)*h(4)*unl(pos4(i)) + h(6)*h(5)*unl(pos5(i)) + h(6)*h(6)*unl(pos6(i))
343 . + h(6)*h(7)*unl(pos7(i)) + h(6)*h(8)*unl(pos8(i)))
344c
345 IF (nodadt == 0) THEN
346 a6 = a6 + vol(i) * xi * (h(6)*h(1)*vnl(pos1(i)) + h(6)*h(2)*vnl(pos2(i)) + h(6)*h(3)*vnl(pos3(i))
347 . + h(6)*h(4)*vnl(pos4(i)) + h(6)*h(5)*vnl(pos5(i)) + h(6)*h(6)*vnl(pos6(i))
348 . + h(6)*h(7)*vnl(pos7(i)) + h(6)*h(8)*vnl(pos8(i)))
349 ELSE
350 a6 = a6 + vol(i) * xi * (h(6)*h(1)*minmasscal*vnl(pos1(i)) +
351 . h(6)*h(2)*minmasscal*vnl(pos2(i)) +
352 . h(6)*h(3)*minmasscal*vnl(pos3(i)) +
353 . h(6)*h(4)*minmasscal*vnl(pos4(i)) +
354 . h(6)*h(5)*minmasscal*vnl(pos5(i)) +
355 . h(6)*h(6)*minmasscal*vnl(pos6(i)) +
356 . h(6)*h(7)*minmasscal*vnl(pos7(i)) +
357 . h(6)*h(8)*minmasscal*vnl(pos8(i)))
358 ENDIF
359c
360 b6 = l2 * vol(i) * ( btb16(i)*unl(pos1(i)) + btb26(i)*unl(pos2(i))
361 . + btb36(i)*unl(pos3(i)) + btb46(i)*unl(pos4(i)) + btb56(i)*unl(pos5(i))
362 . + btb66(i)*unl(pos6(i)) + btb67(i)*unl(pos7(i)) + btb68(i)*unl(pos8(i)))
363c
364 c6 = vol(i) * h(6) * var_reg(i)
365c
366 a7 = vol(i) * (h(7)*h(1)*unl(pos1(i)) + h(7)*h(2)*unl(pos2(i)) + h(7)*h(3)*unl(pos3(i))
367 . + h(7)*h(4)*unl(pos4(i)) + h(7)*h(5)*unl(pos5(i)) + h(7)*h(6)*unl(pos6(i))
368 . + h(7)*h(7)*unl(pos7(i)) + h(7)*h(8)*unl(pos8(i)))
369c
370 IF (nodadt == 0) THEN
371 a7 = a7 + vol(i) * xi * (h(7)*h(1)*vnl(pos1(i)) + h(7)*h(2)*vnl(pos2(i)) + h(7)*h(3)*vnl(pos3(i))
372 . + h(7)*h(4)*vnl(pos4(i)) + h(7)*h(5)*vnl(pos5(i)) + h(7)*h(6)*vnl(pos6(i))
373 . + h(7)*h(7)*vnl(pos7(i)) + h(7)*h(8)*vnl(pos8(i)))
374 ELSE
375 a7 = a7 + vol(i) * xi * (h(7)*h(1)*minmasscal*vnl(pos1(i)) +
376 . h(7)*h(2)*minmasscal*vnl(pos2(i)) +
377 . h(7)*h(3)*minmasscal*vnl(pos3(i)) +
378 . h(7)*h(4)*minmasscal*vnl(pos4(i)) +
379 . h(7)*h(5)*minmasscal*vnl(pos5(i)) +
380 . h(7)*h(6)*minmasscal*vnl(pos6(i)) +
381 . h(7)*h(7)*minmasscal*vnl(pos7(i)) +
382 . h(7)*h(8)*minmasscal*vnl(pos8(i)))
383 ENDIF
384c
385 b7 = l2 * vol(i) * ( btb17(i)*unl(pos1(i)) + btb27(i)*unl(pos2(i))
386 . + btb37(i)*unl(pos3(i)) + btb47(i)*unl(pos4(i)) + btb57(i)*unl(pos5(i))
387 . + btb67(i)*unl(pos6(i)) + btb77(i)*unl(pos7(i)) + btb78(i)*unl(pos8(i)))
388c
389 c7 = vol(i) * h(7) * var_reg(i)
390c
391 a8 = vol(i) * (h(8)*h(1)*unl(pos1(i)) + h(8)*h(2)*unl(pos2(i)) + h(8)*h(3)*unl(pos3(i))
392 . + h(8)*h(4)*unl(pos4(i)) + h(8)*h(5)*unl(pos5(i)) + h(8)*h(6)*unl(pos6(i))
393 . + h(8)*h(7)*unl(pos7(i)) + h(8)*h(8)*unl(pos8(i)))
394c
395 IF (nodadt == 0) THEN
396 a8 = a8 + vol(i) * xi * (h(8)*h(1)*vnl(pos1(i)) + h(8)*h(2)*vnl(pos2(i)) + h(8)*h(3)*vnl(pos3(i))
397 . + h(8)*h(4)*vnl(pos4(i)) + h(8)*h(5)*vnl(pos5(i)) + h(8)*h(6)*vnl(pos6(i))
398 . + h(8)*h(7)*vnl(pos7(i)) + h(8)*h(8)*vnl(pos8(i)))
399 ELSE
400 a8 = a8 + vol(i) * xi * (h(8)*h(1)*minmasscal*vnl(pos1(i)) +
401 . h(8)*h(2)*minmasscal*vnl(pos2(i)) +
402 . h(8)*h(3)*minmasscal*vnl(pos3(i)) +
403 . h(8)*h(4)*minmasscal*vnl(pos4(i)) +
404 . h(8)*h(5)*minmasscal*vnl(pos5(i)) +
405 . h(8)*h(6)*minmasscal*vnl(pos6(i)) +
406 . h(8)*h(7)*minmasscal*vnl(pos7(i)) +
407 . h(8)*h(8)*minmasscal*vnl(pos8(i)))
408 ENDIF
409c
410 b8 = l2 * vol(i) * ( btb18(i)*unl(pos1(i)) + btb28(i)*unl(pos2(i))
411 . + btb38(i)*unl(pos3(i)) + btb48(i)*unl(pos4(i)) + btb58(i)*unl(pos5(i))
412 . + btb68(i)*unl(pos6(i)) + btb78(i)*unl(pos7(i)) + btb88(i)*unl(pos8(i)))
413c
414 c8 = vol(i) * h(8) * var_reg(i)
415c
416c Fint = Vol * (L*L * BT*B*Unl + NT*N*Unl + Damp*NT*N*Vnl - NT*Vreg )
417c The = length stores in the structure
418c Unl = cumul nodal def
419c vnl = velocity non locale
420c
421 f1(i) = a1 + b1 - c1
422 f2(i) = a2 + b2 - c2
423 f3(i) = a3 + b3 - c3
424 f4(i) = a4 + b4 - c4
425 f5(i) = a5 + b5 - c5
426 f6(i) = a6 + b6 - c6
427 f7(i) = a7 + b7 - c7
428 f8(i) = a8 + b8 - c8
429c
430 ! Computing nodal equivalent stiffness
431 IF (nodadt > 0) THEN
432 sti1(i) = (abs(l2*btb11(i) + h(1)*h(1)) + abs(l2*btb12(i) + h(1)*h(2)) + abs(l2*btb13(i) + h(1)*h(3)) +
433 . abs(l2*btb14(i) + h(1)*h(4)) + abs(l2*btb15(i) + h(1)*h(5)) + abs(l2*btb16(i) + h(1)*h(6)) +
434 . abs(l2*btb17(i) + h(1)*h(7)) + abs(l2*btb18(i) + h(1)*h(8)))*vol(i)
435 sti2(i) = (abs(l2*btb12(i) + h(2)*h(1)) + abs(l2*btb22(i) + h(2)*h(2)) + abs(l2*btb23(i) + h(2)*h(3)) +
436 . abs(l2*btb24(i) + h(2)*h(4)) + abs(l2*btb25(i) + h(2)*h(5)) + abs(l2*btb26(i) + h(2)*h(6)) +
437 . abs(l2*btb27(i) + h(2)*h(7)) + abs(l2*btb28(i) + h(2)*h(8)))*vol(i)
438 sti3(i) = (abs(l2*btb13(i) + h(3)*h(1)) + abs(l2*btb23(i) + h(3)*h(2)) + abs(l2*btb33(i) + h(3)*h(3)) +
439 . abs(l2*btb34(i) + h(3)*h(4)) + abs(l2*btb35(i) + h(3)*h(5)) + abs(l2*btb36(i) + h(3)*h(6)) +
440 . abs(l2*btb37(i) + h(3)*h(7)) + abs(l2*btb38(i) + h(3)*h(8)))*vol(i)
441 sti4(i) = (abs(l2*btb14(i) + h(4)*h(1)) + abs(l2*btb24(i) + h(4)*h(2)) + abs(l2*btb34(i) + h(4)*h(3)) +
442 . abs(l2*btb44(i) + h(4)*h(4)) + abs(l2*btb45(i) + h(4)*h(5)) + abs(l2*btb46(i) + h(4)*h(6)) +
443 . abs(l2*btb47(i) + h(4)*h(7)) + abs(l2*btb48(i) + h(4)*h(8)))*vol(i)
444 sti5(i) = (abs(l2*btb15(i) + h(5)*h(1)) + abs(l2*btb25(i) + h(5)*h(2)) + abs(l2*btb35(i) + h(5)*h(3)) +
445 . abs(l2*btb45(i) + h(5)*h(4)) + abs(l2*btb55(i) + h(5)*h(5)) + abs(l2*btb56(i) + h(5)*h(6)) +
446 . abs(l2*btb57(i) + h(5)*h(7)) + abs(l2*btb58(i) + h(5)*h(8)))*vol(i)
447 sti6(i) = (abs(l2*btb16(i) + h(6)*h(1)) + abs(l2*btb26(i) + h(6)*h(2)) + abs(l2*btb36(i) + h(6)*h(3)) +
448 . abs(l2*btb46(i) + h(6)*h(4)) + abs(l2*btb56(i) + h(6)*h(5)) + abs(l2*btb66(i) + h(6)*h(6)) +
449 . abs(l2*btb67(i) + h(6)*h(7)) + abs(l2*btb68(i) + h(6)*h(8)))*vol(i)
450 sti7(i) = (abs(l2*btb17(i) + h(7)*h(1)) + abs(l2*btb27(i) + h(7)*h(2)) + abs(l2*btb37(i) + h(7)*h(3)) +
451 . abs(l2*btb47(i) + h(7)*h(4)) + abs(l2*btb57(i) + h(7)*h(5)) + abs(l2*btb67(i) + h(7)*h(6)) +
452 . abs(l2*btb77(i) + h(7)*h(7)) + abs(l2*btb78(i) + h(7)*h(8)))*vol(i)
453 sti8(i) = (abs(l2*btb18(i) + h(8)*h(1)) + abs(l2*btb28(i) + h(8)*h(2)) + abs(l2*btb38(i) + h(8)*h(3)) +
454 . abs(l2*btb48(i) + h(8)*h(4)) + abs(l2*btb58(i) + h(8)*h(5)) + abs(l2*btb68(i) + h(8)*h(6)) +
455 . abs(l2*btb78(i) + h(8)*h(7)) + abs(l2*btb88(i) + h(8)*h(8)))*vol(i)
456 ENDIF
457c
458 ! If the element is broken
459 ELSE
460c
461 ! Initial element characteristic length
462 lc(i) = ((wi/eight)*vol0(i))**third
463c
464 IF (nodadt > 0) THEN
465 ! Non-local absorbing forces
466 f1(i) = sqrt(mass(pos1(i))/mass0(pos1(i)))*h(1)*zeta*sspnl*half*
467 . (vnl(pos1(i))+vnl0(pos1(i)))*(three/four)*(lc(i)**2)
468 f2(i) = sqrt(mass(pos2(i))/mass0(pos2(i)))*h(2)*zeta*sspnl*half*
469 . (vnl(pos2(i))+vnl0(pos2(i)))*(three/four)*(lc(i)**2)
470 f3(i) = sqrt(mass(pos3(i))/mass0(pos3(i)))*h(3)*zeta*sspnl*half*
471 . (vnl(pos3(i))+vnl0(pos3(i)))*(three/four)*(lc(i)**2)
472 f4(i) = sqrt(mass(pos4(i))/mass0(pos4(i)))*h(4)*zeta*sspnl*half*
473 . (vnl(pos4(i))+vnl0(pos4(i)))*(three/four)*(lc(i)**2)
474 f5(i) = sqrt(mass(pos5(i))/mass0(pos5(i)))*h(5)*zeta*sspnl*half*
475 . (vnl(pos5(i))+vnl0(pos5(i)))*(three/four)*(lc(i)**2)
476 f6(i) = sqrt(mass(pos6(i))/mass0(pos6(i)))*h(6)*zeta*sspnl*half*
477 . (vnl(pos6(i))+vnl0(pos6(i)))*(three/four)*(lc(i)**2)
478 f7(i) = sqrt(mass(pos7(i))/mass0(pos7(i)))*h(7)*zeta*sspnl*half*
479 . (vnl(pos7(i))+vnl0(pos7(i)))*(three/four)*(lc(i)**2)
480 f8(i) = sqrt(mass(pos8(i))/mass0(pos8(i)))*h(8)*zeta*sspnl*half*
481 . (vnl(pos8(i))+vnl0(pos8(i)))*(three/four)*(lc(i)**2)
482 ! Computing nodal equivalent stiffness
483 sti1(i) = em20
484 sti2(i) = em20
485 sti3(i) = em20
486 sti4(i) = em20
487 sti5(i) = em20
488 sti6(i) = em20
489 sti7(i) = em20
490 sti8(i) = em20
491 ELSE
492 ! Non-local absorbing forces
493 f1(i) = h(1)*zeta*sspnl*half*(vnl(pos1(i))+vnl0(pos1(i)))*(three/four)*(lc(i)**2)
494 f2(i) = h(2)*zeta*sspnl*half*(vnl(pos2(i))+vnl0(pos2(i)))*(three/four)*(lc(i)**2)
495 f3(i) = h(3)*zeta*sspnl*half*(vnl(pos3(i))+vnl0(pos3(i)))*(three/four)*(lc(i)**2)
496 f4(i) = h(4)*zeta*sspnl*half*(vnl(pos4(i))+vnl0(pos4(i)))*(three/four)*(lc(i)**2)
497 f5(i) = h(5)*zeta*sspnl*half*(vnl(pos5(i))+vnl0(pos5(i)))*(three/four)*(lc(i)**2)
498 f6(i) = h(6)*zeta*sspnl*half*(vnl(pos6(i))+vnl0(pos6(i)))*(three/four)*(lc(i)**2)
499 f7(i) = h(7)*zeta*sspnl*half*(vnl(pos7(i))+vnl0(pos7(i)))*(three/four)*(lc(i)**2)
500 f8(i) = h(8)*zeta*sspnl*half*(vnl(pos8(i))+vnl0(pos8(i)))*(three/four)*(lc(i)**2)
501 ENDIF
502 ENDIF
503 ENDDO
504c-----------------------------------------------------------------------
505c Assemblage
506c-----------------------------------------------------------------------
507 IF (iparit == 0) THEN
508c
509 fnl => nloc_dmg%FNL(1:l_nloc,itask+1)
510 IF (nodadt > 0) stifnl => nloc_dmg%STIFNL(1:l_nloc,itask+1) ! Non-local equivalent nodal stiffness
511 DO i=1,nel
512 ! Assembling the forces in the classic way
513 fnl(pos1(i)) = fnl(pos1(i)) - f1(i)
514 fnl(pos2(i)) = fnl(pos2(i)) - f2(i)
515 fnl(pos3(i)) = fnl(pos3(i)) - f3(i)
516 fnl(pos4(i)) = fnl(pos4(i)) - f4(i)
517 fnl(pos5(i)) = fnl(pos5(i)) - f5(i)
518 fnl(pos6(i)) = fnl(pos6(i)) - f6(i)
519 fnl(pos7(i)) = fnl(pos7(i)) - f7(i)
520 fnl(pos8(i)) = fnl(pos8(i)) - f8(i)
521 IF (nodadt > 0) THEN
522 ! Spectral radius of stiffness matrix
523 maxstif = max(sti1(i),sti2(i),sti3(i),sti4(i),sti5(i),sti6(i),sti7(i),sti8(i))
524 ! Computing nodal stiffness
525 stifnl(pos1(i)) = stifnl(pos1(i)) + maxstif
526 stifnl(pos2(i)) = stifnl(pos2(i)) + maxstif
527 stifnl(pos3(i)) = stifnl(pos3(i)) + maxstif
528 stifnl(pos4(i)) = stifnl(pos4(i)) + maxstif
529 stifnl(pos5(i)) = stifnl(pos5(i)) + maxstif
530 stifnl(pos6(i)) = stifnl(pos6(i)) + maxstif
531 stifnl(pos7(i)) = stifnl(pos7(i)) + maxstif
532 stifnl(pos8(i)) = stifnl(pos8(i)) + maxstif
533 ENDIF
534 ENDDO
535c
536 ! If PARITH/ON
537 ELSE
538c
539 DO i=1,nel
540 ii = i + nft
541c
542 ! Spectral radius of stiffness matrix
543 IF (nodadt > 0) THEN
544 maxstif = max(sti1(i),sti2(i),sti3(i),sti4(i),sti5(i),sti6(i),sti7(i),sti8(i))
545 ENDIF
546c
547 k = nloc_dmg%IADS(1,ii)
548 IF (ip == 1) THEN
549 nloc_dmg%FSKY(k,1) = -f1(i)
550 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
551 ELSE
552 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f1(i)
553 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
554 ENDIF
555c
556 k = nloc_dmg%IADS(2,ii)
557 IF (ip == 1) THEN
558 nloc_dmg%FSKY(k,1) = -f2(i)
559 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
560 ELSE
561 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f2(i)
562 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
563 ENDIF
564c
565 k = nloc_dmg%IADS(3,ii)
566 IF (ip == 1) THEN
567 nloc_dmg%FSKY(k,1) = -f3(i)
568 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
569 ELSE
570 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f3(i)
571 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
572 ENDIF
573c
574 k = nloc_dmg%IADS(4,ii)
575 IF (ip == 1) THEN
576 nloc_dmg%FSKY(k,1) = -f4(i)
577 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
578 ELSE
579 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f4(i)
580 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
581 ENDIF
582c
583 k = nloc_dmg%IADS(5,ii)
584 IF (ip == 1) THEN
585 nloc_dmg%FSKY(k,1) = -f5(i)
586 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
587 ELSE
588 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f5(i)
589 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
590 ENDIF
591c
592 k = nloc_dmg%IADS(6,ii)
593 IF (ip == 1) THEN
594 nloc_dmg%FSKY(k,1) = -f6(i)
595 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
596 ELSE
597 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f6(i)
598 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
599 ENDIF
600c
601 k = nloc_dmg%IADS(7,ii)
602 IF (ip == 1) THEN
603 nloc_dmg%FSKY(k,1) = -f7(i)
604 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
605 ELSE
606 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f7(i)
607 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
608 ENDIF
609c
610 k = nloc_dmg%IADS(8,ii)
611 IF (ip == 1) THEN
612 nloc_dmg%FSKY(k,1) = -f8(i)
613 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = maxstif
614 ELSE
615 nloc_dmg%FSKY(k,1) = nloc_dmg%FSKY(k,1) - f8(i)
616 IF (nodadt > 0) nloc_dmg%STSKY(k,1) = nloc_dmg%STSKY(k,1) + maxstif
617 ENDIF
618c
619 ENDDO
620 ENDIF
621c
622 !-----------------------------------------------------------------------
623 ! Computing non-local timestep
624 !-----------------------------------------------------------------------
625 IF (nodadt == 0) THEN
626 DO i = 1,nel
627 ! If the element is not broken, normal computation
628 IF (offg(i)/=zero) THEN
629 ! Non-local critical time-step in the plane
630 dtnl = (two*(min(vol(i)**third,le_max))*sqrt(three*zeta))/
631 . sqrt(twelve*l2 + (min(vol(i)**third,le_max))**2)
632 ! Retaining the minimal value
633 dt2t = min(dt2t,dtfac1(1)*cdamp*dtnl)
634 ENDIF
635 ENDDO
636 ENDIF
637c-----------
638 ! Deallocation of tables
639 IF (ALLOCATED(btb11)) DEALLOCATE(btb11)
640 IF (ALLOCATED(btb12)) DEALLOCATE(btb12)
641 IF (ALLOCATED(btb13)) DEALLOCATE(btb13)
642 IF (ALLOCATED(btb14)) DEALLOCATE(btb14)
643 IF (ALLOCATED(btb15)) DEALLOCATE(btb15)
644 IF (ALLOCATED(btb16)) DEALLOCATE(btb16)
645 IF (ALLOCATED(btb17)) DEALLOCATE(btb17)
646 IF (ALLOCATED(btb18)) DEALLOCATE(btb18)
647 IF (ALLOCATED(btb22)) DEALLOCATE(btb22)
648 IF (ALLOCATED(btb23)) DEALLOCATE(btb23)
649 IF (ALLOCATED(btb24)) DEALLOCATE(btb24)
650 IF (ALLOCATED(btb25)) DEALLOCATE(btb25)
651 IF (ALLOCATED(btb26)) DEALLOCATE(btb26)
652 IF (ALLOCATED(btb27)) DEALLOCATE(btb27)
653 IF (ALLOCATED(btb28)) DEALLOCATE(btb28)
654 IF (ALLOCATED(btb33)) DEALLOCATE(btb33)
655 IF (ALLOCATED(btb34)) DEALLOCATE(btb34)
656 IF (ALLOCATED(btb35)) DEALLOCATE(btb35)
657 IF (ALLOCATED(btb36)) DEALLOCATE(btb36)
658 IF (ALLOCATED(btb37)) DEALLOCATE(btb37)
659 IF (ALLOCATED(btb38)) DEALLOCATE(btb38)
660 IF (ALLOCATED(btb44)) DEALLOCATE(btb44)
661 IF (ALLOCATED(btb45)) DEALLOCATE(btb45)
662 IF (ALLOCATED(btb46)) DEALLOCATE(btb46)
663 IF (ALLOCATED(btb47)) DEALLOCATE(btb47)
664 IF (ALLOCATED(btb48)) DEALLOCATE(btb48)
665 IF (ALLOCATED(btb55)) DEALLOCATE(btb55)
666 IF (ALLOCATED(btb56)) DEALLOCATE(btb56)
667 IF (ALLOCATED(btb57)) DEALLOCATE(btb57)
668 IF (ALLOCATED(btb58)) DEALLOCATE(btb58)
669 IF (ALLOCATED(btb66)) DEALLOCATE(btb66)
670 IF (ALLOCATED(btb67)) DEALLOCATE(btb67)
671 IF (ALLOCATED(btb68)) DEALLOCATE(btb68)
672 IF (ALLOCATED(btb77)) DEALLOCATE(btb77)
673 IF (ALLOCATED(btb78)) DEALLOCATE(btb78)
674 IF (ALLOCATED(btb88)) DEALLOCATE(btb88)
675 IF (ALLOCATED(pos1)) DEALLOCATE(pos1)
676 IF (ALLOCATED(pos2)) DEALLOCATE(pos2)
677 IF (ALLOCATED(pos3)) DEALLOCATE(pos3)
678 IF (ALLOCATED(pos4)) DEALLOCATE(pos4)
679 IF (ALLOCATED(pos5)) DEALLOCATE(pos5)
680 IF (ALLOCATED(pos6)) DEALLOCATE(pos6)
681 IF (ALLOCATED(pos7)) DEALLOCATE(pos7)
682 IF (ALLOCATED(pos8)) DEALLOCATE(pos8)
683 IF (ALLOCATED(sti1)) DEALLOCATE(sti1)
684 IF (ALLOCATED(sti2)) DEALLOCATE(sti2)
685 IF (ALLOCATED(sti3)) DEALLOCATE(sti3)
686 IF (ALLOCATED(sti4)) DEALLOCATE(sti4)
687 IF (ALLOCATED(sti5)) DEALLOCATE(sti5)
688 IF (ALLOCATED(sti6)) DEALLOCATE(sti6)
689 IF (ALLOCATED(sti7)) DEALLOCATE(sti7)
690 IF (ALLOCATED(sti8)) DEALLOCATE(sti8)
691c
#define my_real
Definition cppsort.cpp:32
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21