OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i25buc_vox1.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!|| i25buc_vox1 ../starter/source/interfaces/inter3d1/i25buc_vox1.F
25!||--- called by ------------------------------------------------------
26!|| inint3 ../starter/source/interfaces/inter3d1/inint3.F
27!||--- calls -----------------------------------------------------
28!|| ancmsg ../starter/source/output/message/message.F
29!|| i25trivox1 ../starter/source/interfaces/inter3d1/i25trivox1.F
30!|| insol25 ../starter/source/interfaces/inter3d1/i25buc_vox1.F
31!||--- uses -----------------------------------------------------
32!|| message_mod ../starter/share/message_module/message_mod.F
33!|| tri7box ../starter/share/modules1/tri7box.F
34!||====================================================================
35 SUBROUTINE i25buc_vox1(
36 1 X ,IRECT,NSV ,BUMULT,
37 2 NMN ,NRTM ,NSN ,INTBUF_TAB ,
38 3 GAP ,I_STOK ,
39 4 DIST ,TZINF,MAXBOX ,MINBOX,MSR ,
40 5 STF ,STFN ,IDDLEVEL,
41 6 GAP_S,GAP_M ,IGAP ,GAPMIN ,
42 7 GAPMAX,INACTI,GAP_S_L,GAP_M_L,
43 8 MARGE ,ID ,TITR ,NBINFLG,MBINFLG,
44 9 ILEV ,MSEGTYP,GAP_N,BGAPSMX,
45 A IPARTS,KNOD2ELS,NOD2ELS,
46 B KREMNODE,REMNODE,
47 C IXS, IXS10, IXS16, IXS20,ICODE,ISKEW ,
48 D DRAD ,DGAPLOAD,NRTMT,FLAG_REMOVED_NODE,
49 E IELEM_M,nin,npari,ipari )
50C-----------------------------------------------
51C M o d u l e s
52C-----------------------------------------------
53 USE message_mod
54 USE tri7box
56 USE intbufdef_mod
57C-----------------------------------------------
58C I m p l i c i t T y p e s
59C-----------------------------------------------
60#include "implicit_f.inc"
61C-----------------------------------------------
62C C o m m o n B l o c k s
63C-----------------------------------------------
64#include "units_c.inc"
65#include "com04_c.inc"
66#include "scr06_c.inc"
67C-----------------------------------------------
68C D u m m y A r g u m e n t s
69C-----------------------------------------------
70 INTEGER NMN, NRTM, NSN, I_STOK,IGAP,
71 . INACTI,
72 . IPARTS(*), KNOD2ELS(*), NOD2ELS(*)
73 INTEGER IRECT(4,*),NSV(*),ICODE(*),ISKEW(*)
74 INTEGER MSR(*),IDDLEVEL
75 INTEGER NBINFLG(*),MBINFLG(*),ILEV,MSEGTYP(*)
76 INTEGER KREMNODE(*),REMNODE(*)
77 INTEGER IXS(*), IXS10(*), IXS16(*), IXS20(*)
78 LOGICAL, INTENT(in) :: FLAG_REMOVED_NODE !< flag to remove some S node from the list of candidates
79 INTEGER , INTENT(IN) :: IELEM_M(2,NRTM)
80 integer, intent(in) :: nin !< interface index
81 integer, intent(in) :: npari !< first dim of ipari
82 integer, dimension(npari), intent(inout) :: ipari !< interface data
83C REAL
85 . stf(*),stfn(*),x(3,*),gap_s(*),gap_m(*),
86 . dist,bumult,gap,tzinf,maxbox,minbox,gapmin,gapmax,
87 . gap_s_l(*),gap_m_l(*),marge,gap_n(4,*)
89 . bgapsmx
90 my_real , INTENT(IN) :: drad, dgapload
91 INTEGER ID
92 CHARACTER(LEN=NCHARTITLE) :: TITR
93
94 INTEGER , INTENT(IN) :: NRTMT
95 type(intbuf_struct_), intent(inout) :: intbuf_tab !< interface data
96C-----------------------------------------------
97C L o c a l V a r i a b l e s
98C-----------------------------------------------
99 INTEGER I, J, L, N1, N2, N3, N4,N_SOL, ESHIFT
100 INTEGER IBID, NELS
101 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
102C REAL
103 my_real
104 . dx1,dy1,dz1,
105 . dx3,dy3,dz3,
106 . dx4,dy4,dz4,
107 . dx6,dy6,dz6,
108 . dd1,dd2,dd3,dd4,dd,dd0,xmin,ymin,zmin,
109 . xmax,ymax,zmax,tzinf0,minbox0,maxbox0,
110 . bid,tzinf_st,marge_st,dd_st,d_max,pensol,d_moy,
111 . xyzm(6),bminma(6),aaa,ledgmax
112 my_real,
113 . DIMENSION(:), ALLOCATABLE :: edge_l2,edge_l2_tmp
114 INTEGER TAGP(NPART),IAD,N,IE,IL,IP
115 INTEGER, DIMENSION(:), ALLOCATABLE :: NPARTNS,IELEM,LPARTNS
116 INTEGER NBX,NBY,NBZ
117 INTEGER, DIMENSION(:), ALLOCATABLE :: IS_LARGE_NODE,LARGE_NODE, TAGNOD,LOCAL_NEXT_NOD
118 INTEGER :: NB_LARGE_NODES
119 INTEGER (KIND=8) :: NBX8,NBY8,NBZ8,RES8,LVOXEL8
120 INTEGER (KIND=8) :: IONE,IHUNDRED !< Integer constants in INTEGER 8 for comparisions
121C-----------------------------------------------
122 IONE=1
123 ihundred=100
124 gapmax=ep30
125 gapmin=zero
126
127C-------use temporarily GAP_N(1,*)=V/A
128C Not used :: MVOISN(1,*)-> MTYPE(solid),MVOISN(2,*)-> E_id
129C
130C
131C 1-CALCUL TAILLE DES ZONES INFLUENCES
132C DD EST LA LONGEUR MOYENNE ELEMENT
133C DD_ST EST LA LONGEUR MAX ELEMENT
134 dd=zero
135 dd_st=zero
136 pensol=ep30
137 d_moy = zero
138 n_sol = 0
139 DO 10 l=1,nrtm
140C CONNECIVITES ELEMENT
141 n1=irect(1,l)
142 n2=irect(2,l)
143 n3=irect(3,l)
144 n4=irect(4,l)
145C LONGUEUR COTE 1
146 dx1=(x(1,n1)-x(1,n2))
147 dy1=(x(2,n1)-x(2,n2))
148 dz1=(x(3,n1)-x(3,n2))
149 dd1=sqrt(dx1**2+dy1**2+dz1**2)
150C LONGUEUR COTE 2
151 dx3=(x(1,n1)-x(1,n4))
152 dy3=(x(2,n1)-x(2,n4))
153 dz3=(x(3,n1)-x(3,n4))
154 dd2=sqrt(dx3**2+dy3**2+dz3**2)
155C LONGUEUR COTE 3
156 dx4=(x(1,n3)-x(1,n2))
157 dy4=(x(2,n3)-x(2,n2))
158 dz4=(x(3,n3)-x(3,n2))
159 dd3=sqrt(dx4**2+dy4**2+dz4**2)
160C LONGUEUR COTE 4
161 dx6=(x(1,n4)-x(1,n3))
162 dy6=(x(2,n4)-x(2,n3))
163 dz6=(x(3,n4)-x(3,n3))
164 dd4=sqrt(dx6**2+dy6**2+dz6**2)
165 dd=dd+ (dd1+dd2+dd3+dd4)
166C-------only for solid--- and coating shell
167 IF (msegtyp(l)==0.OR.msegtyp(l)>nrtmt) THEN
168 d_max=max(dd1,dd2,dd3,dd4)
169 d_max=min(d_max,gap_n(1,l))
170C--------correction of too huge GAP_N(1,L) w/ irregular mesh
171 gap_n(1,l)=d_max
172 dd_st=max(dd_st,d_max)
173 n_sol = n_sol + 1
174 d_moy = d_moy + d_max
175 END IF
176C IF (MSEGTYP(L)==0) DD_ST=MAX(DD_ST,DD1,DD2,DD3,DD4)
177 10 CONTINUE
178C TAILLE ZINF = .1*TAILLE MOYENNE ELEMENT DE CHAQUE COTE
179C TAILLE BUCKET MIN = TZINF * BUMULT
180 dd0=dd/nrtm/four
181 dd=dd0
182C TZINF = BUMULT*DD
183 marge = bumult*dd
184C calcul de TZINF en fct de la marge et non le contraire
185 tzinf = marge + max(gap+dgapload,drad)
186C MARGE_ST : marge independante du BUMULT en input pour trouver les memes pene initiales (cas delete ou chgt coord.)
187 marge_st = bmul0*dd
188
189C Si IDDLEVEL = 0, alors 1er passage ici, avec la MARGE de l'ENGINE pour trouver les memes candidats
190 IF(iddlevel==0) marge_st = marge
191 tzinf_st = marge_st + max(gap+dgapload,drad)
192 bid = four_over_5*dd
193 IF (inacti/=7.AND.tzinf>bid) THEN
194 ibid = nint(tzinf/dd0)
195 ibid =(2*ibid+4)*ibid*2
196 ENDIF
197C MAXBOX= ZEP9*(DD - TZINF)
198C DD + 2*TZINF : taille element augmentee de zone influence
199 maxbox= half*(dd + 2*tzinf)
200 minbox= half*maxbox
201 tzinf0 = tzinf
202 minbox0 = minbox
203 maxbox0 = maxbox
204C MIS A ZERO POUR FAIRE SEARCH COMPLET CYCLE 0 ENGINE
205 dist = zero
206
207C--------------------------------
208C CALCUL DES BORNES DU DOMAINE
209C--------------------------------
210 bminma(1)=-ep30
211 bminma(2)=-ep30
212 bminma(3)=-ep30
213 bminma(4)= ep30
214 bminma(5)= ep30
215 bminma(6)= ep30
216C
217 xmin=ep30
218 xmax=-ep30
219 ymin=ep30
220 ymax=-ep30
221 zmin=ep30
222 zmax=-ep30
223C
224 DO 20 i=1,nmn
225 j=msr(i)
226 xmin= min(xmin,x(1,j))
227 ymin= min(ymin,x(2,j))
228 zmin= min(zmin,x(3,j))
229 xmax= max(xmax,x(1,j))
230 ymax= max(ymax,x(2,j))
231 zmax= max(zmax,x(3,j))
232 20 CONTINUE
233C
234 xmin=xmin-tzinf_st
235 ymin=ymin-tzinf_st
236 zmin=zmin-tzinf_st
237 xmax=xmax+tzinf_st
238 ymax=ymax+tzinf_st
239 zmax=zmax+tzinf_st
240C
241 DO 25 i=1,nsn
242 j=nsv(i)
243 xmin= min(xmin,x(1,j))
244 ymin= min(ymin,x(2,j))
245 zmin= min(zmin,x(3,j))
246 xmax= max(xmax,x(1,j))
247 ymax= max(ymax,x(2,j))
248 zmax= max(zmax,x(3,j))
249 25 CONTINUE
250C
251 bminma(1) = max(bminma(1),xmax)
252 bminma(2) = max(bminma(2),ymax)
253 bminma(3) = max(bminma(3),zmax)
254 bminma(4) = min(bminma(4),xmin)
255 bminma(5) = min(bminma(5),ymin)
256 bminma(6) = min(bminma(6),zmin)
257
258 xyzm(1) = bminma(4)
259 xyzm(2) = bminma(5)
260 xyzm(3) = bminma(6)
261 xyzm(4) = bminma(1)
262 xyzm(5) = bminma(2)
263 xyzm(6) = bminma(3)
264
265 aaa = sqrt(nmn /
266 . ((bminma(1)-bminma(4))*(bminma(2)-bminma(5))
267 . +(bminma(2)-bminma(5))*(bminma(3)-bminma(6))
268 . +(bminma(3)-bminma(6))*(bminma(1)-bminma(4))))
269
270 aaa = 0.75*aaa
271
272 nbx = nint(aaa*(bminma(1)-bminma(4)))
273 nby = nint(aaa*(bminma(2)-bminma(5)))
274 nbz = nint(aaa*(bminma(3)-bminma(6)))
275 nbx = max(nbx,1)
276 nby = max(nby,1)
277 nbz = max(nbz,1)
278
279 nbx8=nbx
280 nby8=nby
281 nbz8=nbz
282 res8=(nbx8+2)*(nby8+2)*(nbz8+2)
283 lvoxel8 = lvoxel
284
285 IF(res8 > lvoxel8) THEN
286 aaa = lvoxel
287 aaa = aaa/((nbx8+2)*(nby8+2)*(nbz8+2))
288 aaa = aaa**(third)
289 nbx = int((nbx+2)*aaa)-2
290 nby = int((nby+2)*aaa)-2
291 nbz = int((nbz+2)*aaa)-2
292 nbx = max(nbx,1)
293 nby = max(nby,1)
294 nbz = max(nbz,1)
295 ENDIF
296
297 nbx8=nbx
298 nby8=nby
299 nbz8=nbz
300 res8=(nbx8+2)*(nby8+2)*(nbz8+2)
301
302 IF(res8 > lvoxel8) THEN
303 nbx = min(ihundred,max(nbx8,ione))
304 nby = min(ihundred,max(nby8,ione))
305 nbz = min(ihundred,max(nbz8,ione))
306 END IF
307
308C--------------------------------
309C Initializations case of IDDLEVEL=0 or INACTI/=5 and INACTI/=-1
310C--------------------------------
311 ALLOCATE(npartns(nsn+1),ielem(nrtm),edge_l2(nsn))
312
313 npartns(1:nsn+1) = 0
314 ielem(1:nrtm) = 0
315 edge_l2(1:nsn) = zero
316 ledgmax = zero
317C--------------------------------
318
319
320
321 nb_large_nodes = 0
322 ALLOCATE(large_node(nsn))
323 ALLOCATE(is_large_node(nsn))
324 ALLOCATE(tagnod(numnod))
325 is_large_node(1:nsn) = 0
326 large_node(1:nsn) = 0
327 nb_large_nodes = 0
328 tagnod(1:numnod) = 0
329
330 IF(iddlevel==1)THEN
331 ! Looking for initial penetration (vs solids) if and only if IDDLEVEL=1
332C
333 DO ie=1,nrtm
334 ! looks for solid element supporting the segment (cf initial penetrations vs solids)
335 IF(ielem_m(1,ie)<=numels) THEN
336 ielem(ie)= ielem_m(1,ie)
337 ELSE
338 CALL insol25(irect ,ixs ,ixs10,ixs16,ixs20,
339 . knod2els ,nod2els ,nels ,ie )
340 ielem(ie)=nels
341 ENDIF
342 END DO
343C
344 IF(inacti==5.OR.inacti==-1)THEN
345
346 ALLOCATE(edge_l2_tmp(numnod))
347 edge_l2_tmp(1:numnod)=zero
348 n_sol = 0
349 DO ie=1,nrtm
350 IF(stf(ie)> zero)THEN ! Solids only
351 DO il=1,4
352 n1=irect(il,ie)
353 n2=irect(mod(il,4)+1,ie)
354
355 aaa = (x(1,n2)-x(1,n1))*(x(1,n2)-x(1,n1))
356 . + (x(2,n2)-x(2,n1))*(x(2,n2)-x(2,n1))
357 . + (x(3,n2)-x(3,n1))*(x(3,n2)-x(3,n1))
358 edge_l2_tmp(n1) = max(edge_l2_tmp(n1), aaa )
359 edge_l2_tmp(n2) = max(edge_l2_tmp(n2), aaa )
360 IF (msegtyp(ie)==0.OR.msegtyp(ie)>nrtmt) THEN ! Solids only
361 IF(tagnod(n1)==0) THEN
362 n_sol= n_sol + 1
363 tagnod(n1) = 1
364 ENDIF
365 IF(tagnod(n2)==0) THEN
366 n_sol= n_sol + 1
367 tagnod(n2) = 1
368 ENDIF
369 ENDIF
370 END DO
371 ENDIF
372
373 END DO
374C
375
376 DO i=1,nsn
377 n=nsv(i)
378 IF(stfn(i)/=zero) THEN
379 edge_l2(i) = half*sqrt(edge_l2_tmp(n))
380 IF(tagnod(n)==1) ledgmax=ledgmax+edge_l2(i)!MAX(LEDGMAX,EDGE_L2(I))
381 END IF
382 END DO
383
384 IF(n_sol > 0) ledgmax=half*ledgmax/n_sol
385
386
387 DEALLOCATE(edge_l2_tmp)
388
389 ENDIF
390
391c IF(INACTI==5.OR.INACTI==-1)THEN ! Do it for all Inacti ( warnings initial penetrations)
392 tagp(1:npart) =0
393C
394 DO i=1,nsn
395 n=nsv(i)
396 DO iad=knod2els(n)+1,knod2els(n+1)
397 ie=nod2els(iad)
398 ip=iparts(ie)
399 IF(tagp(ip)==0)THEN
400 npartns(i)=npartns(i)+1
401 tagp(ip) =1
402 END IF
403 END DO
404 DO iad=knod2els(n)+1,knod2els(n+1)
405 ie=nod2els(iad)
406 ip=iparts(ie)
407 tagp(ip) =0 ! Reset
408 END DO
409 END DO
410C
411 DO i=1,nsn
412 npartns(i+1) = npartns(i+1) + npartns(i)
413 END DO
414 DO i=nsn,1,-1
415 npartns(i+1) = npartns(i)
416 END DO
417 npartns(1)=0
418C
419 ALLOCATE(lpartns(npartns(nsn+1)))
420C
421 DO i=1,nsn
422 n=nsv(i)
423 DO iad=knod2els(n)+1,knod2els(n+1)
424 ie=nod2els(iad)
425 ip=iparts(ie)
426 IF(tagp(ip)==0)THEN
427 npartns(i)=npartns(i)+1
428 lpartns(npartns(i))=ip
429 tagp(ip) =1
430 END IF
431 END DO
432 DO iad=knod2els(n)+1,knod2els(n+1)
433 ie=nod2els(iad)
434 ip=iparts(ie)
435 tagp(ip) =0 ! Reset
436 END DO
437 END DO
438C
439 DO i=nsn,1,-1
440 npartns(i+1) = npartns(i)
441 END DO
442 npartns(1)=0
443C
444c ELSE
445c ALLOCATE(LPARTNS(0))
446c END IF
447 ELSE
448 ALLOCATE(lpartns(0))
449 END IF
450C
451
452C initialisation complete de VOXEL
453 DO i=inivoxel,(nbx+2)*(nby+2)*(nbz+2)
454 voxel1(i)=0
455 ENDDO
456 inivoxel = max(inivoxel,(nbx+2)*(nby+2)*(nbz+2)+1)
457C
458 eshift=0
459C
460 i_stok = 0
461
462 ALLOCATE(local_next_nod(nsn))
463 ALLOCATE(iix(nsn))
464 ALLOCATE(iiy(nsn))
465 ALLOCATE(iiz(nsn))
466C
467!$OMP PARALLEL
468 CALL i25trivox1(
469 1 nsn ,irect ,x ,
470 2 stfn ,xyzm ,nsv ,i_stok ,
471 3 eshift ,bgapsmx ,
472 4 voxel1 ,nbx ,nby ,nbz ,nrtm ,
473 5 gap_s ,gap_m ,marge_st,
474 6 nbinflg ,mbinflg ,ilev ,msegtyp ,
475 7 igap ,gap_s_l ,gap_m_l ,edge_l2 ,ledgmax ,
476 8 kremnode,remnode,
477 9 iparts ,npartns ,lpartns ,ielem ,icode ,
478 a iskew ,drad, is_large_node, large_node, nb_large_nodes,
479 b dgapload,nrtmt,flag_removed_node,
480 c ielem_m,local_next_nod,iix,iiy,iiz,
481 d intbuf_tab,ipari,nin)
482!$OMP END PARALLEL
483 DEALLOCATE(local_next_nod)
484 DEALLOCATE(iix)
485 DEALLOCATE(iiy)
486 DEALLOCATE(iiz)
487
488 DEALLOCATE(edge_l2,npartns,ielem,lpartns)
489 DEALLOCATE(is_large_node,large_node,tagnod)
490C
491 IF(nsn/=0)THEN
492 WRITE(iout,*)' POSSIBLE IMPACT NUMBER, NSN:',i_stok,nsn
493C
494 ELSE
495 CALL ancmsg(msgid=552,
496 . msgtype=msgwarning,
497 . anmode=aninfo_blind_2,
498 . i1=id,
499 . c1=titr)
500 ENDIF
501C
502 RETURN
503 END
504
505!||====================================================================
506!|| insol25 ../starter/source/interfaces/inter3d1/i25buc_vox1.F
507!||--- called by ------------------------------------------------------
508!|| i25buc_vox1 ../starter/source/interfaces/inter3d1/i25buc_vox1.F
509!||====================================================================
510 SUBROUTINE insol25(IRECT ,IXS ,IXS10,IXS16,IXS20,
511 . KNOD2ELS,NOD2ELS,NEL ,I )
512C-----------------------------------------------
513C I m p l i c i t T y p e s
514C-----------------------------------------------
515#include "implicit_f.inc"
516C-----------------------------------------------
517C C o m m o n B l o c k s
518C-----------------------------------------------
519#include "com04_c.inc"
520C-----------------------------------------------
521C D u m m y A r g u m e n t s
522C-----------------------------------------------
523 INTEGER NEL, I
524 INTEGER IRECT(4,*), IXS(NIXS,*), KNOD2ELS(*), NOD2ELS(*),
525 . ixs10(6,*), ixs16(8,*), ixs20(12,*)
526C-----------------------------------------------
527C L o c a l V a r i a b l e s
528C-----------------------------------------------
529 INTEGER N, JJ, II, K, IC, IAD,
530 . NUSER, NUSERM
531C-----------------------------------------------
532C E x t e r n a l F u n c t i o n s
533C-----------------------------------------------
534C
535 NEL=0
536
537 if(numels==0) RETURN
538
539 IF(irect(1,i)>numnod) RETURN
540
541 ic=0
542 nuserm = -1
543 DO 230 iad=knod2els(irect(1,i))+1,knod2els(irect(1,i)+1)
544 n = nod2els(iad)
545 IF(n <= numels8)THEN
546 DO 210 jj=1,4
547 ii=irect(jj,i)
548 DO k=1,8
549 IF(ixs(k+1,n)==ii) GOTO 210
550 ENDDO
551 GOTO 230
552 210 CONTINUE
553 ic=ic+1
554 nuser = ixs(11,n)
555 IF (nuser>nuserm) THEN
556 nel = n
557 nuserm = nuser
558 ENDIF
559 ELSEIF(n <= numels8+numels10)THEN
560 DO 220 jj=1,4
561 ii=irect(jj,i)
562 DO k=1,8
563 IF(ixs(k+1,n)==ii) GOTO 220
564 ENDDO
565 DO k=1,6
566 IF(ixs10(k,n-numels8)==ii) GOTO 220
567 ENDDO
568 GOTO 230
569 220 CONTINUE
570 ic=ic+1
571 nuser = ixs(11,n)
572 IF (nuser>nuserm) THEN
573 nel = n
574 nuserm = nuser
575 ENDIF
576 ELSEIF(n <= numels8+numels10+numels20)THEN
577 DO 222 jj=1,4
578 ii=irect(jj,i)
579 DO k=1,8
580 IF(ixs(k+1,n)==ii) GOTO 222
581 ENDDO
582 DO k=1,12
583 IF(ixs20(k,n-numels8-numels10)==ii) GOTO 222
584 ENDDO
585 GOTO 230
586 222 CONTINUE
587 ic=ic+1
588 nuser = ixs(11,n)
589 IF (nuser>nuserm) THEN
590 nel = n
591 nuserm = nuser
592 ENDIF
593 ELSEIF(n <= numels8+numels10+numels20+numels16)THEN
594 DO 224 jj=1,4
595 ii=irect(jj,i)
596 DO k=1,8
597 IF(ixs(k+1,n)==ii) GOTO 224
598 ENDDO
599 DO k=1,8
600 IF(ixs16(k,n-numels8-numels10-numels20)==ii) GOTO 224
601 ENDDO
602 GOTO 230
603 224 CONTINUE
604 ic=ic+1
605 nuser = ixs(11,n)
606 IF (nuser>nuserm) THEN
607 nel = n
608 nuserm = nuser
609 ENDIF
610 ELSE
611 GOTO 230
612 END IF
613 230 CONTINUE
614 RETURN
615 END
#define my_real
Definition cppsort.cpp:32
if(complex_arithmetic) id
subroutine i25buc_vox1(x, irect, nsv, bumult, nmn, nrtm, nsn, intbuf_tab, gap, i_stok, dist, tzinf, maxbox, minbox, msr, stf, stfn, iddlevel, gap_s, gap_m, igap, gapmin, gapmax, inacti, gap_s_l, gap_m_l, marge, id, titr, nbinflg, mbinflg, ilev, msegtyp, gap_n, bgapsmx, iparts, knod2els, nod2els, kremnode, remnode, ixs, ixs10, ixs16, ixs20, icode, iskew, drad, dgapload, nrtmt, flag_removed_node, ielem_m, nin, npari, ipari)
Definition i25buc_vox1.F:50
subroutine insol25(irect, ixs, ixs10, ixs16, ixs20, knod2els, nod2els, nel, i)
subroutine i25trivox1(nsn, irect, x, stfn, xyzm, nsv, ii_stok, eshift, bgapsmx, voxel, nbx, nby, nbz, nrtm, gap_s, gap_m, marge, nbinflg, mbinflg, ilev, msegtyp, igap, gap_s_l, gap_m_l, edge_l2, ledgmax, kremnode, remnode, iparts, npartns, lpartns, ielem, icode, iskew, drad, is_large_node, large_node, nb_large_nodes, dgapload, nrtmt, flag_removed_node, ielem_m, local_next_nod, iix, iiy, iiz, intbuf_tab, ipari, nin)
Definition i25trivox1.F:46
subroutine ymax(idn, fac, npc, pld, stiffmin, stiffmax, stiffini, stiffavg)
Definition law100_upd.F:272
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
integer, parameter nchartitle
integer, dimension(lvoxel) voxel1
Definition tri7box.F:53
integer inivoxel
Definition tri7box.F:53
integer lvoxel
Definition tri7box.F:51
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
Definition message.F:889