OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
stat_c_straf.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!|| stat_c_straf ../engine/source/output/sta/stat_c_straf.F
25!||--- called by ------------------------------------------------------
26!|| genstat ../engine/source/output/sta/genstat.F
27!||--- calls -----------------------------------------------------
28!|| spmd_rgather9_dp ../engine/source/mpi/interfaces/spmd_outp.F
29!|| spmd_stat_pgather ../engine/source/mpi/output/spmd_stat.F
30!|| strs_txt50 ../engine/source/output/sta/sta_txt.F
31!|| tab_strs_txt50 ../engine/source/output/sta/sta_txt.F
32!||--- uses -----------------------------------------------------
33!|| elbufdef_mod ../common_source/modules/mat_elem/elbufdef_mod.F90
34!|| my_alloc_mod ../common_source/tools/memory/my_alloc.F90
35!||====================================================================
36 SUBROUTINE stat_c_straf(ELBUF_TAB,IPARG ,IPM ,IGEO ,IXC ,
37 2 IXTG ,WA,WAP0 ,IPARTC, IPARTTG,
38 3 IPART_STATE,STAT_INDXC,STAT_INDXTG,THKE,SIZP0)
39C-----------------------------------------------
40C M o d u l e s
41C-----------------------------------------------
42 USE elbufdef_mod
43 USE my_alloc_mod
44C-----------------------------------------------
45C I m p l i c i t T y p e s
46C-----------------------------------------------
47#include "implicit_f.inc"
48C-----------------------------------------------
49C C o m m o n B l o c k s
50C-----------------------------------------------
51#include "com01_c.inc"
52#include "com04_c.inc"
53#include "param_c.inc"
54#include "units_c.inc"
55#include "scr14_c.inc"
56#include "scr16_c.inc"
57#include "task_c.inc"
58C-----------------------------------------------
59C D u m m y A r g u m e n t s
60C-----------------------------------------------
61 INTEGER SIZLOC,SIZP0
62 INTEGER IXC(NIXC,*),IXTG(NIXTG,*),
63 . iparg(nparg,*),ipm(npropmi,*),igeo(npropgi,*),
64 . ipartc(*), iparttg(*), ipart_state(*),
65 . stat_indxc(*), stat_indxtg(*)
66 my_real
67 . thke(*)
68 TYPE (elbuf_struct_), DIMENSION(NGROUP), TARGET :: elbuf_tab
69 double precision WA(*),WAP0(*)
70C-----------------------------------------------
71C L o c a l V a r i a b l e s
72C-----------------------------------------------
73 INTEGER I,J,K,N,II,JJ,LEN, IOFF, NG, NEL, NFT, ITY, LFT, NPT,
74 . LLT, MLW, ISTRAIN,ID, IPRT0, IPRT,NPG,IPG,IE,NPTR,NPTS,G_STRA,
75 . ithk,kk(8)
76 INTEGER,DIMENSION(:),ALLOCATABLE :: PTWA
77 INTEGER,DIMENSION(:),ALLOCATABLE :: PTWA_P0
78 double precision
79 . thk, em, eb, h1, h2, h3
80 CHARACTER*100 DELIMIT,LINE
81 TYPE(G_BUFEL_) ,POINTER :: GBUF
82 TYPE(l_bufel_) ,POINTER :: LBUF
83 TYPE(buf_lay_) ,POINTER :: BUFLY
84 my_real,
85 . DIMENSION(:),POINTER :: strain
86C-----------------------------------------------
87 DATA delimit(1:60)
88 ./'#---1----|----2----|----3----|----4----|----5----|----6----|'/
89 DATA delimit(61:100)
90 ./'----7----|----8----|----9----|----10---|'/
91C-----------------------------------------------
92C 4-NODE SHELLS
93C-----------------------------------------------
94 CALL my_alloc(ptwa,max(stat_numelc ,stat_numeltg))
95 ALLOCATE(ptwa_p0(0:max(1,stat_numelc_g,stat_numeltg_g)))
96C-----------------------------------------------
97 jj = 0
98 IF(stat_numelc==0) GOTO 200
99
100 ie=0
101 DO ng=1,ngroup
102 ity =iparg(5,ng)
103 IF (ity == 3) THEN
104 gbuf => elbuf_tab(ng)%GBUF
105 mlw =iparg(1,ng)
106 nel =iparg(2,ng)
107 nft =iparg(3,ng)
108 npt = iparg(6,ng)
109 ithk =iparg(28,ng)
110 nptr = elbuf_tab(ng)%NPTR
111 npts = elbuf_tab(ng)%NPTS
112 npg = nptr*npts
113 lft=1
114 llt=nel
115 g_stra = gbuf%G_STRA
116!
117 DO j=1,8 ! length max of GBUF%G_STRA = 8
118 kk(j) = nel*(j-1)
119 ENDDO
120!
121c--------------------
122 DO i=lft,llt
123 n = i + nft
124
125 iprt=ipartc(n)
126 IF(ipart_state(iprt)==0)cycle
127
128 jj = jj + 1
129 IF (mlw /= 0 .AND. mlw /= 13) THEN
130 wa(jj) = gbuf%OFF(i)
131 ELSE
132 wa(jj) = zero
133 ENDIF
134 jj = jj + 1
135 wa(jj) = iprt
136 jj = jj + 1
137 wa(jj) = ixc(nixc,n)
138 jj = jj + 1
139 wa(jj) = npt
140 jj = jj + 1
141 wa(jj) = npg
142 jj = jj + 1
143 IF (mlw /= 0 .AND. mlw /= 13) THEN
144 IF (ithk >0 ) THEN
145 wa(jj) = gbuf%THK(i)
146 ELSE
147 wa(jj) = thke(n)
148 END IF
149 ELSE
150 wa(jj) = zero
151 ENDIF
152c Strain in Gauss points
153 IF (mlw == 0 .or. mlw == 13) THEN
154 DO ipg=1,npg
155 DO j=1,g_stra
156 jj = jj + 1
157 wa(jj)=zero
158 END DO
159 END DO
160 ELSEIF (g_stra /= 0) THEN
161 IF (npg > 1) THEN
162 strain => gbuf%STRPG
163 ELSE
164 strain => gbuf%STRA
165 ENDIF
166 ii = g_stra*(i-1)
167 DO ipg=1,npg
168 k = (ipg-1)*nel*g_stra
169 DO j=1,g_stra
170 jj = jj + 1
171 wa(jj) = strain(kk(j)+i+k)
172 END DO
173 END DO
174 END IF
175
176 ie=ie+1
177C pointeur de fin de zone dans WA
178 ptwa(ie)=jj
179c
180 ENDDO ! I=LFT,LLT
181 END IF ! ITY==3
182 ENDDO ! NG=1,NGROUP
183
184 200 CONTINUE
185
186 IF(nspmd == 1)THEN
187 ptwa_p0(0)=0
188 DO n=1,stat_numelc
189 ptwa_p0(n)=ptwa(n)
190 END DO
191 len=jj
192 DO j=1,len
193 wap0(j)=wa(j)
194 END DO
195 ELSE
196C construit les pointeurs dans le tableau global WAP0
197 CALL spmd_stat_pgather(ptwa,stat_numelc,ptwa_p0,stat_numelc_g)
198 len = 0
199 CALL spmd_rgather9_dp(wa,jj,wap0,sizp0,len)
200 END IF
201
202 IF(ispmd==0.AND.len>0) THEN
203
204 iprt0=0
205 DO n=1,stat_numelc_g
206
207C retrouve le nieme elt dans l'ordre d'id croissant
208 k=stat_indxc(n)
209C retrouve l'adresse dans WAP0
210 j=ptwa_p0(k-1)
211
212 ioff = nint(wap0(j + 1))
213 IF(ioff >= 1)THEN
214 iprt = nint(wap0(j + 2))
215 IF(iprt /= iprt0)THEN
216 IF (izipstrs == 0) THEN
217 WRITE(iugeo,'(A)') delimit
218 WRITE(iugeo,'(A)')'/INISHE/STRA_F'
219 WRITE(iugeo,'(A)')
220 .'#------------------------ REPEAT --------------------------'
221 WRITE(iugeo,'(A)')
222 . '# SHELLID NPT NPG THK'
223 WRITE(iugeo,'(A/A/A)')
224 .'# REPEAT I=1,NPG :',
225 .'# E1, E2, E12, E23, E31,',
226 .'# K1, K2, K12'
227 WRITE(iugeo,'(A)')
228 .'#---------------------- END REPEAT ------------------------'
229 WRITE(iugeo,'(a)') DELIMIT
230 ELSE
231 WRITE(LINE,'(a)') DELIMIT
232 CALL STRS_TXT50(LINE,100)
233 WRITE(LINE,'(a)')'/inishe/stra_f'
234 CALL STRS_TXT50(LINE,100)
235 WRITE(LINE,'(a)')
236 .'#------------------------ REPEAT --------------------------'
237 CALL strs_txt50(line,100)
238 WRITE(line,'(A)')
239 . '# SHELLID NPT NPG THK'
240 CALL strs_txt50(line,100)
241 WRITE(line,'(A)')'# REPEAT I=1,NPG :'
242 CALL strs_txt50(line,100)
243 WRITE(line,'(A)')'# E1, E2, E12, E23, E31,'
244 CALL strs_txt50(line,100)
245 WRITE(line,'(A)')'# K1, K2, K12'
246 CALL strs_txt50(line,100)
247 WRITE(line,'(a)')
248 .'#---------------------- END REPEAT ------------------------'
249 CALL strs_txt50(line,100)
250 WRITE(line,'(A)') delimit
251 CALL strs_txt50(line,100)
252 ENDIF
253 iprt0=iprt
254 END IF
255 id = nint(wap0(j + 3))
256 npt = nint(wap0(j + 4))
257 npg = nint(wap0(j + 5))
258 thk = wap0(j + 6)
259 j = j + 6
260 IF (izipstrs == 0) THEN
261 WRITE(iugeo,'(3I10,1PE20.13)')id,npt,npg,thk
262 ELSE
263 WRITE(line,'(3I10,1PE20.13)')id,npt,npg,thk
264 CALL strs_txt50(line,100)
265 ENDIF
266
267 DO ipg=1,npg
268 IF (izipstrs == 0) THEN
269 WRITE(iugeo,'(1P5E20.13)')(wap0(j + k),k=1,5)
270 WRITE(iugeo,'(1P3E20.13)')(wap0(j + k),k=6,8)
271 ELSE
272 CALL tab_strs_txt50(wap0(1),5,j,sizp0,5)
273 CALL tab_strs_txt50(wap0(6),3,j,sizp0,3)
274 ENDIF
275 END DO
276 END IF
277 ENDDO
278 ENDIF
279
280C-----------------------------------------------
281C 3-NODE SHELLS
282C-----------------------------------------------
283 jj = 0
284 IF (stat_numeltg==0) GOTO 300
285 ie=0
286
287 DO ng=1,ngroup
288 ity =iparg(5,ng)
289 IF (ity == 7) THEN
290 gbuf => elbuf_tab(ng)%GBUF
291 g_stra = gbuf%G_STRA
292 mlw =iparg(1,ng)
293 nel =iparg(2,ng)
294 nft =iparg(3,ng)
295 npt = iparg(6,ng)
296 ithk = iparg(28,ng)
297 nptr = elbuf_tab(ng)%NPTR
298 npts = elbuf_tab(ng)%NPTS
299 npg = nptr*npts
300 lft=1
301 llt=nel
302!
303 DO j=1,8 ! length max of GBUF%G_STRA = 8
304 kk(j) = nel*(j-1)
305 ENDDO
306!
307c--------------------
308 DO i=lft,llt
309 n = i + nft
310
311 iprt=iparttg(n)
312 IF(ipart_state(iprt)==0)cycle
313
314
315 jj = jj + 1
316 IF (mlw /= 0 .AND. mlw /= 13) THEN
317 wa(jj) = gbuf%OFF(i)
318 ELSE
319 wa(jj) = zero
320 ENDIF
321 jj = jj + 1
322 wa(jj) = iprt
323 jj = jj + 1
324 wa(jj) = ixtg(nixtg,n)
325 jj = jj + 1
326 wa(jj) = npt
327 jj = jj + 1
328 wa(jj) = npg
329 jj = jj + 1
330 IF (mlw /= 0 .AND. mlw /= 13) THEN
331 IF (ithk >0 ) THEN
332 wa(jj) = gbuf%THK(i)
333 ELSE
334 wa(jj) = thke(n+numelc)
335 END IF
336 ELSE
337 wa(jj) = zero
338 ENDIF
339
340c Strain in Gauss points
341 IF (mlw == 0 .or. mlw == 13) THEN
342 DO ipg=1,npg
343 DO j=1,g_stra
344 jj = jj + 1
345 wa(jj) = zero
346 END DO
347 END DO
348 ELSEIF (g_stra > 0) THEN
349 IF (npg > 1) THEN
350 strain => gbuf%STRPG
351 ELSE
352 strain => gbuf%STRA
353 ENDIF
354 ii = g_stra*(i-1)
355 DO ipg=1,npg
356 k = (ipg-1)*nel*g_stra
357 DO j=1,g_stra
358 jj = jj + 1
359 wa(jj) = strain(kk(j)+i+k)
360 END DO
361 END DO
362 END IF ! ISTRAIN /=0
363
364 ie=ie+1
365C pointeur de fin de zone
366 ptwa(ie)=jj
367c
368 ENDDO ! I=LFT,LLT
369 END IF ! ITY==3
370 ENDDO ! NG=1,NGROUP
371
372 300 CONTINUE
373
374 IF(nspmd == 1)THEN
375 len=jj
376 DO j=1,len
377 wap0(j)=wa(j)
378 END DO
379 ptwa_p0(0)=0
380 DO n=1,stat_numeltg
381 ptwa_p0(n)=ptwa(n)
382 END DO
383 ELSE
384C construit les pointeurs dans le tableau global WAP0
385 CALL spmd_stat_pgather(ptwa,stat_numeltg,ptwa_p0,stat_numeltg_g)
386 len = 0
387 CALL spmd_rgather9_dp(wa,jj,wap0,sizp0,len)
388 END IF
389
390 IF(ispmd==0.AND.len>0) THEN
391
392 iprt0=0
393 DO n=1,stat_numeltg_g
394
395C retrouve le nieme elt dans l'ordre d'id croissant
396 k=stat_indxtg(n)
397C retrouve l'adresse dans WAP0
398 j=ptwa_p0(k-1)
399
400 ioff = nint(wap0(j + 1))
401 IF(ioff >= 1)THEN
402 iprt = nint(wap0(j + 2))
403 IF(iprt /= iprt0)THEN
404 IF (izipstrs == 0) THEN
405 WRITE(iugeo,'(A)') delimit
406 WRITE(iugeo,'(A)')'/INISH3/STRA_F'
407 WRITE(iugeo,'(A)')
408 .'#------------------------ REPEAT --------------------------'
409 WRITE(iugeo,'(A)')
410 . '# SH3NID NPT NPG THK'
411 WRITE(iugeo,'(A/A/A)')
412 .'# REPEAT I=1,NPG :',
413 .'# E1, E2, E12, E23, E31,',
414 .'# K1, K2, K12'
415 WRITE(iugeo,'(A)')
416 .'#---------------------- END REPEAT ------------------------'
417 WRITE(iugeo,'(A)') delimit
418 ELSE
419 WRITE(line,'(A)') delimit
420 CALL strs_txt50(line,100)
421 WRITE(line,'(A)')'/INISH3/STRA_F'
422 CALL strs_txt50(line,100)
423 WRITE(line,'(A)')
424 .'#------------------------ REPEAT --------------------------'
425 CALL strs_txt50(line,100)
426 WRITE(line,'(A)')
427 . '# SH3NID NPT NPG THK'
428 CALL strs_txt50(line,100)
429 WRITE(line,'(A)')'# REPEAT I=1,NPG :'
430 CALL strs_txt50(line,100)
431 WRITE(line,'(A)')'# E1, E2, E12, E23, E31,'
432 CALL strs_txt50(line,100)
433 WRITE(line,'(A)')'# K1, K2, K12'
434 CALL strs_txt50(line,100)
435 WRITE(line,'(A)')
436 .'#---------------------- END REPEAT ------------------------'
437 CALL strs_txt50(line,100)
438 WRITE(line,'(A)') delimit
439 CALL strs_txt50(line,100)
440 END IF
441 iprt0=iprt
442 END IF
443 id = nint(wap0(j + 3))
444 npt = nint(wap0(j + 4))
445 npg = nint(wap0(j + 5))
446 thk = wap0(j + 6)
447 j = j + 6
448 IF (izipstrs == 0) THEN
449 WRITE(iugeo,'(3I10,1PE20.13)')id,npt,npg,thk
450 ELSE
451 WRITE(line,'(3I10,1PE20.13)')id,npt,npg,thk
452 CALL strs_txt50(line,100)
453 ENDIF
454 DO ipg=1,npg
455 IF (izipstrs == 0) THEN
456 WRITE(iugeo,'(1P5E20.13)')(wap0(j + k),k=1,5)
457 WRITE(iugeo,'(1P3E20.13)')(wap0(j + k),k=6,8)
458 ELSE
459 CALL tab_strs_txt50(wap0(1),5,j,sizp0,5)
460 CALL tab_strs_txt50(wap0(6),3,j,sizp0,3)
461 ENDIF
462 END DO
463 END IF
464
465 ENDDO
466 ENDIF
467c-----------
468 DEALLOCATE(ptwa)
469 DEALLOCATE(ptwa_p0)
470c-----------
471 RETURN
472 END
#define my_real
Definition cppsort.cpp:32
#define max(a, b)
Definition macros.h:21
subroutine spmd_rgather9_dp(v, len, vp0, lenp0, iad)
Definition spmd_outp.F:1015
subroutine spmd_stat_pgather(ptv, ptlen, ptv_p0, ptlen_p0)
Definition spmd_stat.F:53
subroutine strs_txt50(text, length)
Definition sta_txt.F:87
subroutine tab_strs_txt50(wap0, cpt, j, sizp0, nbpline)
Definition sta_txt.F:127
subroutine stat_c_straf(elbuf_tab, iparg, ipm, igeo, ixc, ixtg, wa, wap0, ipartc, iparttg, ipart_state, stat_indxc, stat_indxtg, thke, sizp0)