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

Go to the source code of this file.

Functions/Subroutines

subroutine c3fint_reg (nloc_dmg, var_reg, thk, nel, off, area, nc1, nc2, nc3, px1, py1, py2, bufnl, imat, nddl, itask, dt2t, le, thk0, area0, nft)

Function/Subroutine Documentation

◆ c3fint_reg()

subroutine c3fint_reg ( type(nlocal_str_), target nloc_dmg,
intent(inout) var_reg,
intent(in) thk,
integer nel,
intent(in) off,
intent(in) area,
integer, dimension(nel) nc1,
integer, dimension(nel) nc2,
integer, dimension(nel) nc3,
intent(in) px1,
intent(in) py1,
intent(in) py2,
type(buf_nloc_), target bufnl,
integer imat,
integer nddl,
integer itask,
intent(inout) dt2t,
intent(in) le,
intent(in) thk0,
intent(in) area0,
integer, intent(in) nft )

Definition at line 31 of file c3fint_reg.F.

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,px1,py1,py2,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
91 ! Safety coefficient for non-local stability vs mechanical stability
92 ! (it has been slightly increased vs nloc_dmg_init.F)
93 my_real, PARAMETER :: csta = 40.0d0
94 ! Coefficient for non-local stability to take into account damping
95 my_real, PARAMETER :: cdamp = 0.7d0
96c-----------------------------------------------------------------------
97c VAR_REG : variable a regulariser (local cumulated plastic strain)
98c NTVAR = NT * VAR_REG
99C=======================================================================
100 ! Recovering non-local parameters
101 l2 = nloc_dmg%LEN(imat)**2 ! Non-local internal length ** 2
102 xi = nloc_dmg%DAMP(imat) ! Non-local damping parameter
103 zeta = nloc_dmg%DENS(imat) ! Non-local density
104 sspnl = nloc_dmg%SSPNL(imat) ! Non-local sound speed
105 l_nloc = nloc_dmg%L_NLOC ! Length of non-local tables
106 le_max = nloc_dmg%LE_MAX(imat) ! Maximal length of convergence
107 ! Allocation of elementary forces vectors
108 ALLOCATE(f1(nel,nddl),f2(nel,nddl),f3(nel,nddl))
109 ! Only for nodal timestep
110 IF (nodadt > 0) THEN
111 ! Non-local nodal stiffness
112 ALLOCATE(sti1(nel,nddl),sti2(nel,nddl),sti3(nel,nddl))
113 ! Non-local mass
114 mass => nloc_dmg%MASS(1:l_nloc)
115 ! Initial non-local mass
116 mass0 => nloc_dmg%MASS0(1:l_nloc)
117 ENDIF
118 ALLOCATE(btb11(nel),btb12(nel),btb13(nel),btb22(nel),
119 . btb23(nel),btb33(nel),vol(nel),pos1(nel),
120 . pos2(nel),pos3(nel))
121 ! Recovering non-local data
122 vnl => nloc_dmg%VNL(1:l_nloc) ! Non-local variable velocities
123 vnl0 => nloc_dmg%VNL_OLD(1:l_nloc) ! Non-local variable velocities
124 unl => nloc_dmg%UNL(1:l_nloc) ! Non-local cumulated variable
125 ! Constant for triangle elements
126 ntn = three*three
127c
128 !-----------------------------------------------------------------------
129 ! Computation of the element volume and the BtB matrix product
130 !-----------------------------------------------------------------------
131 ! Loop over elements
132# include "vectorize.inc"
133 DO i=1,nel
134c
135 ! Recovering the nodes of the triangle element
136 n1 = nloc_dmg%IDXI(nc1(i))
137 n2 = nloc_dmg%IDXI(nc2(i))
138 n3 = nloc_dmg%IDXI(nc3(i))
139c
140 ! Recovering the positions of the first d.o.fs of each nodes
141 pos1(i) = nloc_dmg%POSI(n1)
142 pos2(i) = nloc_dmg%POSI(n2)
143 pos3(i) = nloc_dmg%POSI(n3)
144c
145 ! Computation of the element volume
146 vol(i) = area(i)*thk(i)
147c
148 ! Computation of the product LEN**2 * BtxB
149 btb11(i) = px1(i)**2 + py1(i)**2
150 btb12(i) = -px1(i)**2 + py1(i)*py2(i)
151 btb13(i) = -py1(i)*(py1(i)+py2(i))
152 btb22(i) = px1(i)**2 + py2(i)**2
153 btb23(i) = -py2(i)*(py1(i)+py2(i))
154 btb33(i) = (py1(i)+py2(i))**2
155c
156 ENDDO
157c
158 !-----------------------------------------------------------------------
159 ! Pre-treatment con-local regularization in the thickness
160 !-----------------------------------------------------------------------
161 ! Only if NDDL > 1
162 IF ((nddl > 1).AND.(l2>zero)) THEN
163c
164 ! Allocation of the velocities predictor
165 IF (nddl > 2) THEN
166 IF (nodadt > 0) THEN
167 ALLOCATE(stifnlth(nel,nddl+1))
168 ALLOCATE(dtn(nel,nddl+1))
169 ENDIF
170 ndnod = nddl+1
171 ELSE
172 IF (nodadt > 0) THEN
173 ALLOCATE(stifnlth(nel,nddl))
174 ALLOCATE(dtn(nel,nddl))
175 ENDIF
176 ndnod = nddl
177 ENDIF
178c
179 ! Pointing the non-local values in the thickness of the corresponding element
180 massth => bufnl%MASSTH(1:nel,1:ndnod)
181 unlth => bufnl%UNLTH(1:nel ,1:ndnod)
182 vnlth => bufnl%VNLTH(1:nel ,1:ndnod)
183 fnlth => bufnl%FNLTH(1:nel ,1:ndnod)
184c
185 DO k = 1,ndnod
186 DO i = 1,nel
187 ! Resetting non-local forces
188 fnlth(i,k) = zero
189 ! Resetting non-local nodal stiffness
190 IF (nodadt > 0) THEN
191 stifnlth(i,k) = em20
192 ENDIF
193 ENDDO
194 ENDDO
195c
196 ! Computation of non-local forces in the shell thickness
197 DO k = 1, nddl
198c
199 ! Computation of shape functions value
200 IF ((nddl==2).AND.(k==2)) THEN
201 nth1 = (z0(k,nddl) - zth(k,nddl)) / (zth(k-1,nddl) - zth(k,nddl))
202 nth2 = (z0(k,nddl) - zth(k-1,nddl)) / (zth(k,nddl) - zth(k-1,nddl))
203 ELSE
204 nth1 = (z0(k,nddl) - zth(k+1,nddl)) / (zth(k,nddl) - zth(k+1,nddl))
205 nth2 = (z0(k,nddl) - zth(k,nddl)) / (zth(k+1,nddl) - zth(k,nddl))
206 ENDIF
207c
208 ! Loop over elements
209 DO i = 1,nel
210 ! Computation of B-matrix values
211 IF ((nddl==2).AND.(k==2)) THEN
212 bth1 = (one/(zth(k-1,nddl) - zth(k,nddl)))*(one/thk(i))
213 bth2 = (one/(zth(k,nddl) - zth(k-1,nddl)))*(one/thk(i))
214 ELSE
215 bth1 = (one/(zth(k,nddl) - zth(k+1,nddl)))*(one/thk(i))
216 bth2 = (one/(zth(k+1,nddl) - zth(k,nddl)))*(one/thk(i))
217 ENDIF
218c
219 ! Computation of the non-local K matrix
220 k1 = l2*(bth1**2) + nth1**2
221 k12 = l2*(bth1*bth2)+ (nth1*nth2)
222 k2 = l2*(bth2**2) + nth2**2
223c
224 ! Computation of the non-local forces
225 IF ((nddl==2).AND.(k==2)) THEN
226 fnlth(i,k-1) = fnlth(i,k-1) + (k1*unlth(i,k-1) + k12*unlth(i,k)
227 . + xi*((nth1**2)*vnlth(i,k-1)
228 . + (nth1*nth2)*vnlth(i,k))
229 . - (nth1*var_reg(i,k)))*vol(i)*wf(k,nddl)
230 fnlth(i,k) = fnlth(i,k) + (k12*unlth(i,k-1) + k2*unlth(i,k)
231 . + xi*(nth1*nth2*vnlth(i,k-1)
232 . + (nth2**2)*vnlth(i,k))
233 . - nth2*var_reg(i,k))*vol(i)*wf(k,nddl)
234 ELSE
235 fnlth(i,k) = fnlth(i,k) + (k1*unlth(i,k) + k12*unlth(i,k+1)
236 . + xi*((nth1**2)*vnlth(i,k)
237 . + (nth1*nth2)*vnlth(i,k+1))
238 . - (nth1*var_reg(i,k)))*vol(i)*wf(k,nddl)
239 fnlth(i,k+1) = fnlth(i,k+1) + (k12*unlth(i,k) + k2*unlth(i,k+1)
240 . + xi*(nth1*nth2*vnlth(i,k)
241 . + (nth2**2)*vnlth(i,k+1))
242 . - nth2*var_reg(i,k))*vol(i)*wf(k,nddl)
243 ENDIF
244c
245 ! Computation of non-local nodal stiffness
246 IF (nodadt > 0) THEN
247 IF ((nddl==2).AND.(k==2)) THEN
248 stifnlth(i,k-1) = stifnlth(i,k-1) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
249 stifnlth(i,k) = stifnlth(i,k) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
250 ELSE
251 stifnlth(i,k) = stifnlth(i,k) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
252 stifnlth(i,k+1) = stifnlth(i,k+1) + max(abs(k1)+abs(k12),abs(k12)+abs(k2))*vol(i)*wf(k,nddl)
253 ENDIF
254 ENDIF
255c
256 ENDDO
257 ENDDO
258c
259 ! Updating non-local mass with /DT/NODA
260 IF (nodadt > 0) THEN
261C
262 ! Initial computation of the nodal timestep
263 dtnod = ep20
264 DO k = 1,ndnod
265 DO i = 1,nel
266 dtn(i,k) = dtfac1(11)*cdamp*sqrt(two * massth(i,k) / max(stifnlth(i,k),em20))
267 dtnod = min(dtn(i,k),dtnod)
268 ENDDO
269 ENDDO
270C
271 ! /DT/NODA/CSTX - Constant timestep with added mass
272 IF ((idtmin(11)==3).OR.(idtmin(11)==4).OR.(idtmin(11)==8)) THEN
273 ! Added mass computation if necessary
274 IF (dtnod < dtmin1(11)) THEN
275 DO k = 1,ndnod
276 DO i = 1,nel
277 IF (dtn(i,k) < dtmin1(11)) THEN
278 dt2p = dtmin1(11)/(dtfac1(11)*cdamp)
279 massth(i,k) = max(massth(i,k),csta*half*stifnlth(i,k)*dt2p*dt2p*onep00001)
280 ENDIF
281 ENDDO
282 ENDDO
283 ENDIF
284 dtnod = dtmin1(11)*(sqrt(csta))
285 ENDIF
286C
287 ! Classical nodal timestep check
288 IF (dtnod < dt2t) THEN
289 dt2t = min(dt2t,dtnod)
290 ENDIF
291 ENDIF
292c
293 DO k = 1,ndnod
294 DO i = 1,nel
295 ! Updating the non-local in-thickness velocities
296 vnlth(i,k) = vnlth(i,k) - (fnlth(i,k)/massth(i,k))*dt12
297 ENDDO
298 ENDDO
299c
300 DO k = 1,ndnod
301 DO i = 1,nel
302 ! Computing the non-local in-thickness cumulated values
303 unlth(i,k) = unlth(i,k) + vnlth(i,k)*dt1
304 ENDDO
305 ENDDO
306c
307 ! Transfert at the integration point
308 DO k = 1, nddl
309 !Computation of shape functions value
310 IF ((nddl==2).AND.(k==2)) THEN
311 nth1 = (z0(k,nddl) - zth(k,nddl))/(zth(k-1,nddl) - zth(k,nddl))
312 nth2 = (z0(k,nddl) - zth(k-1,nddl)) /(zth(k,nddl) - zth(k-1,nddl))
313 ELSE
314 nth1 = (z0(k,nddl) - zth(k+1,nddl))/(zth(k,nddl) - zth(k+1,nddl))
315 nth2 = (z0(k,nddl) - zth(k,nddl)) /(zth(k+1,nddl) - zth(k,nddl))
316 ENDIF
317 ! Loop over elements
318 DO i = 1,nel
319 !Integration points non-local variables
320 IF ((nddl==2).AND.(k==2)) THEN
321 var_reg(i,k) = nth1*unlth(i,k-1) + nth2*unlth(i,k)
322 ELSE
323 var_reg(i,k) = nth1*unlth(i,k) + nth2*unlth(i,k+1)
324 ENDIF
325 ENDDO
326 ENDDO
327 ENDIF
328c
329 !-----------------------------------------------------------------------
330 ! Computation of the elementary non-local forces
331 !-----------------------------------------------------------------------
332 ! Loop over additional degrees of freedom
333 DO k = 1, nddl
334c
335 ! Loop over elements
336# include "vectorize.inc"
337 DO i = 1, nel
338c
339 ! If the element is not broken, normal computation
340 IF (off(i) /= zero) THEN
341 ! Computation of the product LEN**2 * BtxB
342 ! Warning: the derivatives of the shape function does not take into account the volume of the element.
343 ! That is why a factor (1/VOL**2) is added in B1, B2, B3.
344 b1 = (l2 / vol(i)) * wf(k,nddl)*(btb11(i)*unl(pos1(i)+k-1) + btb12(i)*unl(pos2(i)+k-1)
345 . + btb13(i)*unl(pos3(i)+k-1))
346c
347 b2 = (l2 / vol(i)) * wf(k,nddl)*(btb12(i)*unl(pos1(i)+k-1) + btb22(i)*unl(pos2(i)+k-1)
348 . + btb23(i)*unl(pos3(i)+k-1))
349c
350 b3 = (l2 / vol(i)) * wf(k,nddl)*(btb13(i)*unl(pos1(i)+k-1) + btb23(i)*unl(pos2(i)+k-1)
351 . + btb33(i)*unl(pos3(i)+k-1))
352c
353 ! Computing the product NtN*UNL
354 ntn_unl = (unl(pos1(i)+k-1) + unl(pos2(i)+k-1) + unl(pos3(i)+k-1))/ntn
355c
356 ! Computing the product XDAMP*NtN*VNL
357 ntn_vnl = (vnl(pos1(i)+k-1) + vnl(pos2(i)+k-1) + vnl(pos3(i)+k-1))/ntn
358 IF (nodadt > 0) THEN
359 ntn_vnl = min(sqrt(mass(pos1(i)+k-1)/mass0(pos1(i)+k-1)),
360 . sqrt(mass(pos2(i)+k-1)/mass0(pos2(i)+k-1)),
361 . sqrt(mass(pos3(i)+k-1)/mass0(pos3(i)+k-1)))*ntn_vnl
362 ENDIF
363c
364 ! Multiplication by the volume of the element
365 ntn_unl = ntn_unl * vol(i) * wf(k,nddl)
366 ntn_vnl = ntn_vnl * xi * vol(i) * wf(k,nddl)
367c
368 ! Introducing the internal variable to be regularized
369 ntvar = var_reg(i,k)*third*vol(i)*wf(k,nddl)
370c
371 ! Computing the elementary non-local forces
372 a = ntn_unl + ntn_vnl - ntvar
373 f1(i,k) = a + b1
374 f2(i,k) = a + b2
375 f3(i,k) = a + b3
376c
377 ! Computing nodal equivalent stiffness
378 IF (nodadt > 0) THEN
379 sti1(i,k) = wf(k,nddl)*(abs((l2/vol(i))*btb11(i) + one/ntn*vol(i)) +
380 . abs((l2/vol(i))*btb12(i) + one/ntn*vol(i)) +
381 . abs((l2/vol(i))*btb13(i) + one/ntn*vol(i)))
382 sti2(i,k) = wf(k,nddl)*(abs((l2/vol(i))*btb12(i) + one/ntn*vol(i)) +
383 . abs((l2/vol(i))*btb22(i) + one/ntn*vol(i)) +
384 . abs((l2/vol(i))*btb23(i) + one/ntn*vol(i)))
385 sti3(i,k) = wf(k,nddl)*(abs((l2/vol(i))*btb13(i) + one/ntn*vol(i)) +
386 . abs((l2/vol(i))*btb23(i) + one/ntn*vol(i)) +
387 . abs((l2/vol(i))*btb33(i) + one/ntn*vol(i)))
388 ENDIF
389c
390 ! If the element is broken, the non-local wave is absorbed
391 ELSE
392 IF (nodadt > 0) THEN
393 ! Non-local absorbing forces
394 f1(i,k) = wf(k,nddl)*sqrt(mass(pos1(i)+k-1)/mass0(pos1(i)+k-1))*zeta*sspnl*
395 . half*(vnl(pos1(i)+k-1)+vnl0(pos1(i)+k-1))*sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
396 f2(i,k) = wf(k,nddl)*sqrt(mass(pos2(i)+k-1)/mass0(pos2(i)+k-1))*zeta*sspnl*
397 . half*(vnl(pos2(i)+k-1)+vnl0(pos2(i)+k-1))*sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
398 f3(i,k) = wf(k,nddl)*sqrt(mass(pos3(i)+k-1)/mass0(pos3(i)+k-1))*zeta*sspnl*
399 . half*(vnl(pos3(i)+k-1)+vnl0(pos3(i)+k-1))*sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
400 ! Computing nodal equivalent stiffness
401 sti1(i,k) = em20
402 sti2(i,k) = em20
403 sti3(i,k) = em20
404 ELSE
405 ! Non-local absorbing forces
406 f1(i,k) = wf(k,nddl)*zeta*sspnl*half*(vnl(pos1(i)+k-1)+vnl0(pos1(i)+k-1))*
407 . sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
408 f2(i,k) = wf(k,nddl)*zeta*sspnl*half*(vnl(pos2(i)+k-1)+vnl0(pos2(i)+k-1))*
409 . sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
410 f3(i,k) = wf(k,nddl)*zeta*sspnl*half*(vnl(pos3(i)+k-1)+vnl0(pos3(i)+k-1))*
411 . sqrt((four/sqrt(three))*(area0(i)))*thk0(i)
412 ENDIF
413 ENDIF
414 ENDDO
415 ENDDO
416c
417 !-----------------------------------------------------------------------
418 ! Assembling of the non-local forces
419 !-----------------------------------------------------------------------
420c
421 ! If PARITH/OFF
422 IF (iparit == 0) THEN
423 ! Recovering non-local internal forces
424 fnl => nloc_dmg%FNL(1:l_nloc,itask+1) ! Non-local forces
425 IF (nodadt > 0) stifnl => nloc_dmg%STIFNL(1:l_nloc,itask+1) ! Non-local equivalent nodal stiffness
426 ! Loop over elements
427 DO i=1,nel
428 ! Loop over non-local degrees of freedom (do not switch the two loops)
429#include "vectorize.inc"
430 DO k = 1,nddl
431 ! Assembling non-local forces
432 fnl(pos1(i)+k-1) = fnl(pos1(i)+k-1) - f1(i,k)
433 fnl(pos2(i)+k-1) = fnl(pos2(i)+k-1) - f2(i,k)
434 fnl(pos3(i)+k-1) = fnl(pos3(i)+k-1) - f3(i,k)
435 IF (nodadt > 0) THEN
436 ! Spectral radius of stiffness matrix
437 maxstif = max(sti1(i,k),sti2(i,k),sti3(i,k))
438 ! Computing nodal stiffness
439 stifnl(pos1(i)+k-1) = stifnl(pos1(i)+k-1) + maxstif
440 stifnl(pos2(i)+k-1) = stifnl(pos2(i)+k-1) + maxstif
441 stifnl(pos3(i)+k-1) = stifnl(pos3(i)+k-1) + maxstif
442 ENDIF
443 ENDDO
444 ENDDO
445c
446 ! If PARITH/ON
447 ELSE
448 ! Loop over additional d.o.fs
449 DO j = 1,nddl
450c
451 ! Loop over elements
452 DO i=1,nel
453 ii = i + nft
454c
455 ! Spectral radius of stiffness matrix
456 IF (nodadt > 0) THEN
457 maxstif = max(sti1(i,j),sti2(i,j),sti3(i,j))
458 ENDIF
459c
460 k = nloc_dmg%IADTG(1,ii)
461 nloc_dmg%FSKY(k,j) = -f1(i,j)
462 IF (nodadt > 0) nloc_dmg%STSKY(k,j) = maxstif
463c
464 k = nloc_dmg%IADTG(2,ii)
465 nloc_dmg%FSKY(k,j) = -f2(i,j)
466 IF (nodadt > 0) nloc_dmg%STSKY(k,j) = maxstif
467c
468 k = nloc_dmg%IADTG(3,ii)
469 nloc_dmg%FSKY(k,j) = -f3(i,j)
470 IF (nodadt > 0) nloc_dmg%STSKY(k,j) = maxstif
471c
472 ENDDO
473c
474 ENDDO
475 ENDIF
476c
477 !-----------------------------------------------------------------------
478 ! Computing non-local timestep
479 !-----------------------------------------------------------------------
480 IF (nodadt == 0) THEN
481 DO i = 1,nel
482 ! If the element is not broken, normal computation
483 IF (off(i)/=zero) THEN
484 ! Non-local critical time-step in the plane
485 dtnl = (two*(min(le(i),le_max))*sqrt(three*zeta))/
486 . sqrt(twelve*l2 + (min(le(i),le_max))**2)
487 ! Non-local critical time-step in the thickness
488 IF (nddl>1) THEN
489 IF (nddl > 2) THEN
490 dtnl_th = (two*(min(thk(i)/nddl,le_max))*sqrt(three*zeta))/
491 . sqrt(twelve*l2 + (min(thk(i)/nddl,le_max))**2)
492 ELSE
493 dtnl_th = (two*(min(thk(i),le_max))*sqrt(three*zeta))/
494 . sqrt(twelve*l2 + (min(thk(i),le_max))**2)
495 ENDIF
496 ELSE
497 dtnl_th = ep20
498 ENDIF
499 ! Retaining the minimal value
500 dt2t = min(dt2t,dtfac1(1)*cdamp*dtnl_th,dtfac1(1)*cdamp*dtnl)
501 ENDIF
502 ENDDO
503 ENDIF
504c
505 ! Deallocation of tables
506 IF (ALLOCATED(f1)) DEALLOCATE(f1)
507 IF (ALLOCATED(f2)) DEALLOCATE(f2)
508 IF (ALLOCATED(f3)) DEALLOCATE(f3)
509 IF (ALLOCATED(sti1)) DEALLOCATE(sti1)
510 IF (ALLOCATED(sti2)) DEALLOCATE(sti2)
511 IF (ALLOCATED(sti3)) DEALLOCATE(sti3)
512 IF (ALLOCATED(btb11)) DEALLOCATE(btb11)
513 IF (ALLOCATED(btb12)) DEALLOCATE(btb12)
514 IF (ALLOCATED(btb13)) DEALLOCATE(btb13)
515 IF (ALLOCATED(btb22)) DEALLOCATE(btb22)
516 IF (ALLOCATED(btb23)) DEALLOCATE(btb23)
517 IF (ALLOCATED(btb33)) DEALLOCATE(btb33)
518 IF (ALLOCATED(stifnlth)) DEALLOCATE(stifnlth)
519 IF (ALLOCATED(pos1)) DEALLOCATE(pos1)
520 IF (ALLOCATED(pos2)) DEALLOCATE(pos2)
521 IF (ALLOCATED(pos3)) DEALLOCATE(pos3)
522 IF (ALLOCATED(vol)) DEALLOCATE(vol)
523c
#define my_real
Definition cppsort.cpp:32
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