OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
cdk6fint_reg.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23!||====================================================================
24!|| cdk6fint_reg ../engine/source/elements/sh3n/coquedk6/cdk6fint_reg.F
25!||--- called by ------------------------------------------------------
26!|| cdk6forc3 ../engine/source/elements/sh3n/coquedk6/cdk6forc3.F
27!||--- uses -----------------------------------------------------
28!|| elbufdef_mod ../common_source/modules/mat_elem/elbufdef_mod.F90
29!|| nlocal_reg_mod ../common_source/modules/nlocal_reg_mod.F
30!||====================================================================
31 SUBROUTINE cdk6fint_reg(
32 1 NLOC_DMG,VAR_REG, THK, NEL,
33 2 OFF, AREA, NC1, NC2,
34 3 NC3, PX2, PY2, PX3,
35 4 PY3, BUFNL, IMAT, NDDL,
36 5 ITASK, DT2T, LE, THK0,
37 6 AREA0, NFT)
38C-----------------------------------------------
39C M o d u l e s
40C-----------------------------------------------
42 USE elbufdef_mod
43C-----------------------------------------------
44C I m p l i c i t T y p e s
45C-----------------------------------------------
46#include "implicit_f.inc"
47C-----------------------------------------------
48C G l o b a l P a r a m e t e r s
49C-----------------------------------------------
50#include "parit_c.inc"
51#include "com20_c.inc"
52#include "com08_c.inc"
53#include "scr02_c.inc"
54#include "scr18_c.inc"
55C-----------------------------------------------
56C D u m m y A r g u m e n t s
57C-----------------------------------------------
58 INTEGER, INTENT(IN) :: NFT
59 INTEGER :: NEL,IMAT,NDDL,ITASK
60 INTEGER, DIMENSION(NEL) :: NC1,NC2,NC3
61 my_real, DIMENSION(NEL,NDDL), INTENT(INOUT)::
62 . VAR_REG
63 my_real, DIMENSION(NEL), INTENT(IN) ::
64 . area,off,px2,py2,px3,py3,thk,le,thk0,area0
65 my_real, INTENT(INOUT) ::
66 . dt2t
67 TYPE(nlocal_str_), TARGET :: NLOC_DMG
68 TYPE(BUF_NLOC_) , TARGET :: BUFNL
69C-----------------------------------------------
70C L o c a l V a r i a b l e s
71C-----------------------------------------------
72 INTEGER I,K,N1,N2,N3,L_NLOC,II,J,NDNOD
73 my_real
74 . l2,ntn,ntn_unl,ntn_vnl,xi,ntvar,a,
75 . b1,b2,b3,zeta,sspnl,
76 . nth1,nth2,bth1,bth2,k1,k2,k12,
77 . dtnl_th,dtnl,le_max,maxstif,
78 . dtnod,dt2p
79 my_real, DIMENSION(:,:), ALLOCATABLE ::
80 . f1,f2,f3,sti1,sti2,sti3
81 my_real, DIMENSION(:) ,ALLOCATABLE ::
82 . btb11,btb12,btb13,btb22,btb23,btb33,vol
83 INTEGER, DIMENSION(:), ALLOCATABLE ::
84 . POS1,POS2,POS3
85 my_real, POINTER, DIMENSION(:) ::
86 . vnl,fnl,unl,stifnl,mass,mass0,vnl0
87 my_real, POINTER, DIMENSION(:,:) ::
88 . massth,unlth,vnlth,fnlth
89 my_real, DIMENSION(:,:), ALLOCATABLE ::
90 . stifnlth,dtn
91c-----------------------------------------------------------------------
92c VAR_REG : variable a regulariser (local cumulated plastic strain)
93c NTVAR = NT * VAR_REG
94C=======================================================================
95 ! Recovering non-local parameters
96 l2 = nloc_dmg%LEN(imat)**2 ! Non-local internal length ** 2
97 xi = nloc_dmg%DAMP(imat) ! Non-local damping parameter
98 zeta = nloc_dmg%DENS(imat) ! Non-local density
99 sspnl = nloc_dmg%SSPNL(imat) ! Non-local sound speed
100 l_nloc = nloc_dmg%L_NLOC ! Length of non-local tables
101 le_max = nloc_dmg%LE_MAX(imat) ! Maximal length of convergence
102 ! Allocation of elementary forces vectors
103 ALLOCATE(f1(nel,nddl),f2(nel,nddl),f3(nel,nddl))
104 ! Only for nodal timestep
105 IF (nodadt > 0) THEN
106 ! Non-local nodal stiffness
107 ALLOCATE(sti1(nel,nddl),sti2(nel,nddl),sti3(nel,nddl))
108 ! Non-local mass
109 mass => nloc_dmg%MASS(1:l_nloc)
110 ! Initial non-local mass
111 mass0 => nloc_dmg%MASS0(1:l_nloc)
112 ENDIF
113 ALLOCATE(btb11(nel),btb12(nel),btb13(nel),btb22(nel),
114 . btb23(nel),btb33(nel),vol(nel),pos1(nel),
115 . pos2(nel),pos3(nel))
116 ! Recovering non-local data
117 vnl => nloc_dmg%VNL(1:l_nloc) ! Non-local variable velocities
118 vnl0 => nloc_dmg%VNL_OLD(1:l_nloc) ! Non-local variable velocities
119 unl => nloc_dmg%UNL(1:l_nloc) ! Non-local cumulated variable
120 ! Constant for triangle elements
121 ntn = three*three
122c
123 !-----------------------------------------------------------------------
124 ! Computation of the element volume and the BtB matrix product
125 !-----------------------------------------------------------------------
126 ! Loop over elements
127# include "vectorize.inc"
128 DO i=1,nel
129c
130 ! Recovering the nodes of the triangle element
131 n1 = nloc_dmg%IDXI(nc1(i))
132 n2 = nloc_dmg%IDXI(nc2(i))
133 n3 = nloc_dmg%IDXI(nc3(i))
134c
135 ! Recovering the positions of the first d.o.fs of each nodes
136 pos1(i) = nloc_dmg%POSI(n1)
137 pos2(i) = nloc_dmg%POSI(n2)
138 pos3(i) = nloc_dmg%POSI(n3)
139c
140 ! Computation of the element volume
141 vol(i) = area(i)*thk(i)
142c
143 ! Computation of the product BtxB
144 btb11(i) = (-px3(i)-px2(i))**2 + (-py2(i)-py3(i))**2
145 btb12(i) = (-px3(i)-px2(i))*px2(i) + (-py2(i)-py3(i))*py2(i)
146 btb13(i) = (-px3(i)-px2(i))*px3(i) + (-py2(i)-py3(i))*py3(i)
147 btb22(i) = px2(i)**2 + py2(i)**2
148 btb23(i) = px2(i)*px3(i) + py2(i)*py3(i)
149 btb33(i) = px3(i)**2 + py3(i)**2
150c
151 ENDDO
152c
153 !-----------------------------------------------------------------------
154 ! Pre-treatment non-local regularization in the thickness
155 !-----------------------------------------------------------------------
156 ! Only if NDDL > 1
157 IF ((nddl > 1).AND.(l2>zero)) THEN
158c
159 ! Allocation of the velocities predictor
160 IF (nddl > 2) THEN
161 IF (nodadt > 0) THEN
162 ALLOCATE(stifnlth(nel,nddl+1))
163 ALLOCATE(dtn(nel,nddl+1))
164 ENDIF
165 ndnod = nddl+1
166 ELSE
167 IF (nodadt > 0) THEN
168 ALLOCATE(stifnlth(nel,nddl))
169 ALLOCATE(dtn(nel,nddl))
170 ENDIF
171 ndnod = nddl
172 ENDIF
173c
174 ! Pointing the non-local values in the thickness of the corresponding element
175 massth => bufnl%MASSTH(1:nel,1:ndnod)
176 unlth => bufnl%UNLTH(1:nel,1:ndnod)
177 vnlth => bufnl%VNLTH(1:nel,1:ndnod)
178 fnlth => bufnl%FNLTH(1:nel,1:ndnod)
179c
180 DO k = 1,ndnod
181 DO i = 1,nel
182 ! Resetting non-local forces
183 fnlth(i,k) = zero
184 ! Resetting non-local nodal stiffness
185 IF (nodadt > 0) THEN
186 stifnlth(i,k) = em20
187 ENDIF
188 ENDDO
189 ENDDO
190c
191 ! Computation of non-local forces in the shell thickness
192 DO k = 1, nddl
193c
194 ! Computation of shape functions value
195 IF ((nddl==2).AND.(k==2)) THEN
196 nth1 = (z0(k,nddl) - zth(k,nddl)) / (zth(k-1,nddl) - zth(k,nddl))
197 nth2 = (z0(k,nddl) - zth(k-1,nddl)) / (zth(k,nddl) - zth(k-1,nddl))
198 ELSE
199 nth1 = (z0(k,nddl) - zth(k+1,nddl)) / (zth(k,nddl) - zth(k+1,nddl))
200 nth2 = (z0(k,nddl) - zth(k,nddl)) / (zth(k+1,nddl) - zth(k,nddl))
201 ENDIF
202c
203 ! Loop over elements
204 DO i = 1,nel
205 ! Computation of B-matrix values
206 IF ((nddl==2).AND.(k==2)) THEN
207 bth1 = (one/(zth(k-1,nddl) - zth(k,nddl)))*(one/thk(i))
208 bth2 = (one/(zth(k,nddl) - zth(k-1,nddl)))*(one/thk(i))
209 ELSE
210 bth1 = (one/(zth(k,nddl) - zth(k+1,nddl)))*(one/thk(i))
211 bth2 = (one/(zth(k+1,nddl) - zth(k,nddl)))*(one/thk(i))
212 ENDIF
213c
214 ! Computation of the non-local K matrix
215 k1 = l2*(bth1**2) + nth1**2
216 k12 = l2*(bth1*bth2)+ (nth1*nth2)
217 k2 = l2*(bth2**2) + nth2**2
218c
219 ! Computation of the non-local forces
220 IF ((nddl==2).AND.(k==2)) THEN
221 fnlth(i,k-1) = fnlth(i,k-1) + (k1*unlth(i,k-1) + k12*unlth(i,k)
222 . + xi*((nth1**2)*vnlth(i,k-1)
223 . + (nth1*nth2)*vnlth(i,k))
224 . - (nth1*var_reg(i,k)))*vol(i)*wf(k,nddl)
225 fnlth(i,k) = fnlth(i,k) + (k12*unlth(i,k-1) + k2*unlth(i,k)
226 . + xi*(nth1*nth2*vnlth(i,k-1)
227 . + (nth2**2)*vnlth(i,k))
228 . - nth2*var_reg(i,k))*vol(i)*wf(k,nddl)
229 ELSE
230 fnlth(i,k) = fnlth(i,k) + (k1*unlth(i,k) + k12*unlth(i,k+1)
231 . + xi*((nth1**2)*vnlth(i,k)
232 . + (nth1*nth2)*vnlth(i,k+1))
233 . - (nth1*var_reg(i,k)))*vol(i)*wf(k,nddl)
234 fnlth(i,k+1) = fnlth(i,k+1) + (k12*unlth(i,k) + k2*unlth(i,k+1)
235 . + xi*(nth1*nth2*vnlth(i,k)
236 . + (nth2**2)*vnlth(i,k+1))
237 . - nth2*var_reg(i,k))*vol(i)*wf(k,nddl)
238 ENDIF
239c
240 ! Computation of non-local nodal stiffness
241 IF (nodadt > 0) THEN
242 IF ((nddl==2).AND.(k==2)) THEN
243 stifnlth(i,k-1) = stifnlth(i,k-1) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
244 stifnlth(i,k) = stifnlth(i,k) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
245 ELSE
246 stifnlth(i,k) = stifnlth(i,k) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
247 stifnlth(i,k+1) = stifnlth(i,k+1) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
248 ENDIF
249 ENDIF
250c
251 ENDDO
252 ENDDO
253c
254 ! Updating non-local mass with /DT/NODA
255 IF (nodadt > 0) THEN
256C
257 ! Initial computation of the nodal timestep
258 dtnod = ep20
259 DO k = 1,ndnod
260 DO i = 1,nel
261 dtn(i,k) = dtfac1(11)*cdamp*sqrt(two*massth(i,k)/max(stifnlth(i,k),em20))
262 dtnod = min(dtn(i,k),dtnod)
263 ENDDO
264 ENDDO
265C
266 ! /DT/NODA/CSTX - Constant timestep with added mass
267 IF ((idtmin(11)==3).OR.(idtmin(11)==4).OR.(idtmin(11)==8)) THEN
268 ! Added mass computation if necessary
269 IF (dtnod < dtmin1(11)*(sqrt(csta))) THEN
270 DO k = 1,ndnod
271 DO i = 1,nel
272 IF (dtn(i,k) < dtmin1(11)) THEN
273 dt2p = dtmin1(11)/(dtfac1(11)*cdamp)
274 massth(i,k) = max(massth(i,k),csta*half*stifnlth(i,k)*dt2p*dt2p*onep00001)
275 ENDIF
276 ENDDO
277 ENDDO
278 dtnod = dtmin1(11)*(sqrt(csta))
279 ENDIF
280 ENDIF
281C
282 ! Classical nodal timestep check
283 IF (dtnod < dt2t) THEN
284 dt2t = min(dt2t,dtnod)
285 ENDIF
286 ENDIF
287c
288 DO k = 1,ndnod
289 DO i = 1,nel
290 ! Updating the non-local in-thickness velocities
291 vnlth(i,k) = vnlth(i,k) - (fnlth(i,k)/massth(i,k))*dt12
292 ENDDO
293 ENDDO
294c
295 DO k = 1,ndnod
296 DO i = 1,nel
297 ! Computing the non-local in-thickness cumulated values
298 unlth(i,k) = unlth(i,k) + vnlth(i,k)*dt1
299 ENDDO
300 ENDDO
301c
302 ! Transfert at the integration point
303 DO k = 1, nddl
304 !Computation of shape functions value
305 IF ((nddl==2).AND.(k==2)) THEN
306 nth1 = (z0(k,nddl) - zth(k,nddl))/(zth(k-1,nddl) - zth(k,nddl))
307 nth2 = (z0(k,nddl) - zth(k-1,nddl)) /(zth(k,nddl) - zth(k-1,nddl))
308 ELSE
309 nth1 = (z0(k,nddl) - zth(k+1,nddl))/(zth(k,nddl) - zth(k+1,nddl))
310 nth2 = (z0(k,nddl) - zth(k,nddl)) /(zth(k+1,nddl) - zth(k,nddl))
311 ENDIF
312 ! Loop over elements
313 DO i = 1,nel
314 !Integration points non-local variables
315 IF ((nddl==2).AND.(k==2)) THEN
316 var_reg(i,k) = nth1*unlth(i,k-1) + nth2*unlth(i,k)
317 ELSE
318 var_reg(i,k) = nth1*unlth(i,k) + nth2*unlth(i,k+1)
319 ENDIF
320 ENDDO
321 ENDDO
322 ENDIF
323c
324 !-----------------------------------------------------------------------
325 ! Computation of the elementary non-local forces
326 !-----------------------------------------------------------------------
327 ! Loop over additional degrees of freedom
328 DO k = 1, nddl
329c
330 ! Loop over elements
331# include "vectorize.inc"
332 DO i=1,nel
333c
334 ! If the element is not broken, normal computation
335 IF (off(i)/=zero) THEN
336 ! Computation of the product LEN**2 * BtxB
337 b1 = (l2 * vol(i)) * wf(k,nddl)*(btb11(i)*unl(pos1(i)+k-1) + btb12(i)*unl(pos2(i)+k-1)
338 . + btb13(i)*unl(pos3(i)+k-1))
339c
340 b2 = (l2 * vol(i)) * wf(k,nddl)*(btb12(i)*unl(pos1(i)+k-1) + btb22(i)*unl(pos2(i)+k-1)
341 . + btb23(i)*unl(pos3(i)+k-1))
342c
343 b3 = (l2 * vol(i)) * wf(k,nddl)*(btb13(i)*unl(pos1(i)+k-1) + btb23(i)*unl(pos2(i)+k-1)
344 . + btb33(i)*unl(pos3(i)+k-1))
345c
346 ! Computing the product NtN*UNL
347 ntn_unl = (unl(pos1(i)+k-1) + unl(pos2(i)+k-1) + unl(pos3(i)+k-1))/ntn
348c
349 ! Computing the product XDAMP*NtN*VNL
350 ntn_vnl = (vnl(pos1(i)+k-1) + vnl(pos2(i)+k-1) + vnl(pos3(i)+k-1))/ntn
351 IF (nodadt > 0) THEN
352 ntn_vnl = min(sqrt(mass(pos1(i)+k-1)/mass0(pos1(i)+k-1)),
353 . sqrt(mass(pos2(i)+k-1)/mass0(pos2(i)+k-1)),
354 . sqrt(mass(pos3(i)+k-1)/mass0(pos3(i)+k-1)))*ntn_vnl
355 ENDIF
356c
357 ! Multiplication by the volume of the element
358 ntn_unl = ntn_unl * vol(i) * wf(k,nddl)
359 ntn_vnl = ntn_vnl * xi * vol(i) * wf(k,nddl)
360c
361 ! Introducing the internal variable to be regularized
362 ntvar = var_reg(i,k)*third*vol(i)*wf(k,nddl)
363c
364 ! Computing the elementary non-local forces
365 a = ntn_unl + ntn_vnl - ntvar
366 f1(i,k) = a + b1
367 f2(i,k) = a + b2
368 f3(i,k) = a + b3
369c
370 ! Computing nodal equivalent stiffness
371 IF (nodadt > 0) THEN
372 sti1(i,k) = (abs(l2*btb11(i) + one/ntn) + abs(l2*btb12(i) + one/ntn) +
373 . abs(l2*btb13(i) + one/ntn))* vol(i) * wf(k,nddl)
374 sti2(i,k) = (abs(l2*btb12(i) + one/ntn) + abs(l2*btb22(i) + one/ntn) +
375 . abs(l2*btb23(i) + one/ntn))* vol(i) * wf(k,nddl)
376 sti3(i,k) = (abs(l2*btb13(i) + one/ntn) + abs(l2*btb23(i) + one/ntn) +
377 . abs(l2*btb33(i) + one/ntn))* vol(i) * wf(k,nddl)
378 ENDIF
379c
380 ! If the element is broken, the non-local wave is absorbed
381 ELSE
382 IF (nodadt > 0) THEN
383 ! Non-local absorbing forces
384 f1(i,k) = wf(k,nddl)*sqrt(mass(pos1(i)+k-1)/mass0(pos1(i)+k-1))*zeta*sspnl*
385 . half*(vnl(pos1(i)+k-1)+vnl0(pos1(i)+k-1))*sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
386 f2(i,k) = wf(k,nddl)*sqrt(mass(pos2(i)+k-1)/mass0(pos2(i)+k-1))*zeta*sspnl*
387 . half*(vnl(pos2(i)+k-1)+vnl0(pos2(i)+k-1))*sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
388 f3(i,k) = wf(k,nddl)*sqrt(mass(pos3(i)+k-1)/mass0(pos3(i)+k-1))*zeta*sspnl*
389 . half*(vnl(pos3(i)+k-1)+vnl0(pos3(i)+k-1))*sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
390 ! Computing nodal equivalent stiffness
391 sti1(i,k) = em20
392 sti2(i,k) = em20
393 sti3(i,k) = em20
394 ELSE
395 ! Non-local absorbing forces
396 f1(i,k) = wf(k,nddl)*zeta*sspnl*half*(vnl(pos1(i)+k-1)+vnl0(pos1(i)+k-1))*
397 . sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
398 f2(i,k) = wf(k,nddl)*zeta*sspnl*half*(vnl(pos2(i)+k-1)+vnl0(pos2(i)+k-1))*
399 . sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
400 f3(i,k) = wf(k,nddl)*zeta*sspnl*half*(vnl(pos3(i)+k-1)+vnl0(pos3(i)+k-1))*
401 . sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
402 ENDIF
403 ENDIF
404 ENDDO
405 ENDDO
406c
407 !-----------------------------------------------------------------------
408 ! Assembling of the non-local forces
409 !-----------------------------------------------------------------------
410c
411 ! If PARITH/OFF
412 IF (iparit == 0) THEN
413 ! Recovering non-local internal forces
414 fnl => nloc_dmg%FNL(1:l_nloc,itask+1) ! Non-local forces
415 IF (nodadt > 0) stifnl => nloc_dmg%STIFNL(1:l_nloc,itask+1) ! Non-local equivalent nodal stiffnes
416 ! Loop over elements
417 DO i=1,nel
418 ! Loop over non-local degrees of freedom (do not switch the two loops)
419#include "vectorize.inc"
420 DO k = 1,nddl
421 ! Assembling non-local forces
422 fnl(pos1(i)+k-1) = fnl(pos1(i)+k-1) - f1(i,k)
423 fnl(pos2(i)+k-1) = fnl(pos2(i)+k-1) - f2(i,k)
424 fnl(pos3(i)+k-1) = fnl(pos3(i)+k-1) - f3(i,k)
425 IF (nodadt > 0) THEN
426 ! Spectral radius of stiffness matrix
427 maxstif = max(sti1(i,k),sti2(i,k),sti3(i,k))
428 ! Computing nodal stiffness
429 stifnl(pos1(i)+k-1) = stifnl(pos1(i)+k-1) + maxstif
430 stifnl(pos2(i)+k-1) = stifnl(pos2(i)+k-1) + maxstif
431 stifnl(pos3(i)+k-1) = stifnl(pos3(i)+k-1) + maxstif
432 ENDIF
433 ENDDO
434 ENDDO
435c
436 ! If PARITH/ON
437 ELSE
438 ! Loop over additional d.o.fs
439 DO j = 1,nddl
440c
441 ! Loop over elements
442 DO i=1,nel
443 ii = i + nft
444c
445 ! Spectral radius of stiffness matrix
446 IF (nodadt > 0) THEN
447 maxstif = max(sti1(i,j),sti2(i,j),sti3(i,j))
448 ENDIF
449c
450 k = nloc_dmg%IADTG(1,ii)
451 nloc_dmg%FSKY(k,j) = -f1(i,j)
452 IF (nodadt > 0) nloc_dmg%STSKY(k,j) = maxstif
453c
454 k = nloc_dmg%IADTG(2,ii)
455 nloc_dmg%FSKY(k,j) = -f2(i,j)
456 IF (nodadt > 0) nloc_dmg%STSKY(k,j) = maxstif
457c
458 k = nloc_dmg%IADTG(3,ii)
459 nloc_dmg%FSKY(k,j) = -f3(i,j)
460 IF (nodadt > 0) nloc_dmg%STSKY(k,j) = maxstif
461c
462 ENDDO
463c
464 ENDDO
465 ENDIF
466c
467 !-----------------------------------------------------------------------
468 ! Computing non-local timestep
469 !-----------------------------------------------------------------------
470 IF (nodadt == 0) THEN
471 DO i = 1,nel
472 ! If the element is not broken, normal computation
473 IF (off(i)/=zero) THEN
474 ! Non-local critical time-step in the plane
475 dtnl = (two*(min(le(i),le_max))*sqrt(three*zeta))/
476 . sqrt(twelve*l2 + (min(le(i),le_max))**2)
477 ! Non-local critical time-step in the thickness
478 IF (nddl>1) THEN
479 IF (nddl > 2) THEN
480 dtnl_th = (two*(min(thk(i)/nddl,le_max))*sqrt(three*zeta))/
481 . sqrt(twelve*l2 + (min(thk(i)/nddl,le_max))**2)
482 ELSE
483 dtnl_th = (two*(min(thk(i),le_max))*sqrt(three*zeta))/
484 . sqrt(twelve*l2 + (min(thk(i),le_max))**2)
485 ENDIF
486 ELSE
487 dtnl_th = ep20
488 ENDIF
489 ! Retaining the minimal value
490 dt2t = min(dt2t,dtfac1(1)*cdamp*dtnl_th,dtfac1(1)*cdamp*dtnl)
491 ENDIF
492 ENDDO
493 ENDIF
494c
495 ! Deallocation of tables
496 IF (ALLOCATED(f1)) DEALLOCATE(f1)
497 IF (ALLOCATED(f2)) DEALLOCATE(f2)
498 IF (ALLOCATED(f3)) DEALLOCATE(f3)
499 IF (ALLOCATED(sti1)) DEALLOCATE(sti1)
500 IF (ALLOCATED(sti2)) DEALLOCATE(sti2)
501 IF (ALLOCATED(sti3)) DEALLOCATE(sti3)
502 IF (ALLOCATED(btb11)) DEALLOCATE(btb11)
503 IF (ALLOCATED(btb12)) DEALLOCATE(btb12)
504 IF (ALLOCATED(btb13)) DEALLOCATE(btb13)
505 IF (ALLOCATED(btb22)) DEALLOCATE(btb22)
506 IF (ALLOCATED(btb23)) DEALLOCATE(btb23)
507 IF (ALLOCATED(btb33)) DEALLOCATE(btb33)
508 IF (ALLOCATED(stifnlth)) DEALLOCATE(stifnlth)
509 IF (ALLOCATED(pos1)) DEALLOCATE(pos1)
510 IF (ALLOCATED(pos2)) DEALLOCATE(pos2)
511 IF (ALLOCATED(pos3)) DEALLOCATE(pos3)
512 IF (ALLOCATED(vol)) DEALLOCATE(vol)
513 END
subroutine cdk6fint_reg(nloc_dmg, var_reg, thk, nel, off, area, nc1, nc2, nc3, px2, py2, px3, py3, bufnl, imat, nddl, itask, dt2t, le, thk0, area0, nft)
subroutine area(d1, x, x2, y, y2, eint, stif0)
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21