OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
intvo8.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!|| intvo8 ../engine/source/interfaces/inter3d/intvo8.F
25!||--- called by ------------------------------------------------------
26!|| intfop8 ../engine/source/interfaces/interf/intfop8.F
27!||--- calls -----------------------------------------------------
28!|| i8cor3 ../engine/source/interfaces/inter3d/i8cor3.F
29!|| i8cst3 ../engine/source/interfaces/inter3d/i8cst3.F
30!|| i8dis3 ../engine/source/interfaces/inter3d/i8dis3.F
31!|| i8for3 ../engine/source/interfaces/inter3d/i8for3.F
32!|| i8gap3 ../engine/source/interfaces/inter3d/i8gap3.F
33!|| i8loc3 ../engine/source/interfaces/inter3d/i8loc3.F
34!|| i8msr3 ../engine/source/interfaces/inter3d/i8msr3.F
35!|| spmd_i8_commslv ../engine/source/mpi/interfaces/spmd_i8tool.F
36!|| spmd_i8_iloc ../engine/source/mpi/interfaces/spmd_i8tool.F
37!|| spmd_i8_index ../engine/source/mpi/interfaces/spmd_i8tool.F
38!|| spmd_i8_irtl ../engine/source/mpi/interfaces/spmd_i8tool.F
39!|| spmd_i8_updbuf ../engine/source/mpi/interfaces/spmd_i8tool.F
40!||--- uses -----------------------------------------------------
41!|| h3d_mod ../engine/share/modules/h3d_mod.f
42!|| int8_mod ../common_source/modules/interfaces/int8_mod.F90
43!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
44!||====================================================================
45 SUBROUTINE intvo8(IPARI,X ,A ,
46 2 ICODT,FSAV ,V ,MS ,
47 3 FSKYI,ISKY ,FCONT,FNCONT,FTCONT,
48 4 ICONTACT,RCONTACT,
49 5 STIFN,ITAB,INTBUF_TAB, T8, H3D_DATA ,
50 6 NIN , PSKIDS,TAGNCONT,KLOADPINTER,LOADPINTER,
51 7 LOADP_HYD_INTER)
52C-----------------------------------------------
53C M o d u l e s
54C-----------------------------------------------
55 USE intbufdef_mod
56 USE int8_mod
57 USE h3d_mod
58C-----------------------------------------------
59C I m p l i c i t T y p e s
60C-----------------------------------------------
61#include "implicit_f.inc"
62C-----------------------------------------------
63C C o m m o n B l o c k s
64C-----------------------------------------------
65#include "com01_c.inc"
66#include "com04_c.inc"
67#include "com08_c.inc"
68#include "parit_c.inc"
69C-----------------------------------------------
70C D u m m y A r g u m e n t s
71C-----------------------------------------------
72 INTEGER IPARI(*), ICODT(*), ISKY(*),
73 . ICONTACT(*) ,NIN
74 INTEGER ITAB(*),
75 . TAGNCONT(NLOADP_HYD_INTER,NUMNOD),
76 . KLOADPINTER(NINTER+1),LOADPINTER(NINTER*NLOADP_HYD),
77 . LOADP_HYD_INTER(NLOADP_HYD)
78C REAL
80 . x(*), a(3,*), fsav(*), v(3,*), ms(*),
81 . fskyi(lskyi,nfskyi),fcont(3,*), fncont(3,*), ftcont(3,*),
82 . rcontact(*),stifn(*),pskids(*)
83
84 TYPE(intbuf_struct_) INTBUF_TAB
85 TYPE(int8_struct_) T8
86 TYPE(H3D_DATABASE) :: H3D_DATA
87C-----------------------------------------------
88C L o c a l V a r i a b l e s
89C-----------------------------------------------
90 INTEGER I,SIZ
91 INTEGER, DIMENSION(:,:), ALLOCATABLE ::
92 . tab_rmax_uid
93 INTEGER, DIMENSION(:), ALLOCATABLE ::
94 . HAS_MOVED,
95 . INDEX_IN_COMM,
96 . ix1, ix2, ix3,ix4,
97 . nsv_a, !NSV reduced only to active secnd on the proc
98 . nsv2, !Numbering from NSV_A to NSV
99 . irtl_a !IRTL reduced only to active secnd on the proc
100C REAL
101 my_real, DIMENSION(:), ALLOCATABLE :: tab_rmax
102 my_real, DIMENSION(:), ALLOCATABLE :: distance,dist
103 my_real, DIMENSION(:), ALLOCATABLE :: x1,y1,z1,x2,y2,z2,x3,y3,z3
104 my_real, DIMENSION(:), ALLOCATABLE :: x4,y4,z4
105 my_real, DIMENSION(:), ALLOCATABLE :: xi,xi_a,yi,yi_a,zi,zi_a
106 my_real, DIMENSION(:), ALLOCATABLE :: xp,yp,zp,n1,n2,n3
107 my_real, DIMENSION(:), ALLOCATABLE :: ans,ssc,ttc
108 my_real, DIMENSION(:), ALLOCATABLE :: h1,h2,h3,h4
109 my_real, DIMENSION(:), ALLOCATABLE :: face,face_a,stif
110 my_real, DIMENSION(:), ALLOCATABLE :: x0,xx1,xx2,xx3,xx4
111 my_real, DIMENSION(:), ALLOCATABLE :: y0,yy1,yy2,yy3,yy4
112 my_real, DIMENSION(:), ALLOCATABLE :: z0,zz1,zz2,zz3,zz4
113 my_real, DIMENSION(:), ALLOCATABLE :: xi1,xi2,xi3,xi4
114 my_real, DIMENSION(:), ALLOCATABLE :: yi1,yi2,yi3,yi4
115 my_real, DIMENSION(:), ALLOCATABLE :: zi1,zi2,zi3,zi4
116 my_real, DIMENSION(:), ALLOCATABLE :: xn1,xn2,xn3,xn4
117 my_real, DIMENSION(:), ALLOCATABLE :: yn1,yn2,yn3,yn4
118 my_real, DIMENSION(:), ALLOCATABLE :: zn1,zn2,zn3,zn4
119 my_real, DIMENSION(:), ALLOCATABLE :: thk,area,alp
120 my_real, DIMENSION(:), ALLOCATABLE :: distlin
121 INTEGER IBC, NCIMP, NRTM4, NOINT,IGIMP,NINSKID
122 INTEGER MFROT,IBAG,IADM,IFORM,IFT0,IFLINEAR
123C REAL
124 my_real
125 . startt, fric, gap, stopt,
126 . visc,viscf,fnor,depth,fric_last,fnor_last
127 INTEGER :: NMN, NTY, NSN
128 INTEGER :: LFT, LLT, NFT
129C-----------------------------------------------
130
131 nsn = ipari(5)
132 nmn = ipari(6)
133 nty = ipari(7)
134 ibc =ipari(11)
135 ncimp =ipari(13)
136 nrtm4 =ipari(14)
137 noint =ipari(15)
138 mfrot =ipari(30)
139 ibag =ipari(32)
140 iadm =ipari(44)
141 iform=ipari(48)
142 ift0 =ipari(50)
143 iflinear = ipari(49)
144C
145 startt=intbuf_tab%VARIABLES(3)
146 stopt =intbuf_tab%VARIABLES(11)
147 IF(startt>tt) RETURN
148 IF(tt>stopt) RETURN
149C
150 fric =intbuf_tab%VARIABLES(1)
151 gap =intbuf_tab%VARIABLES(2)
152 visc =intbuf_tab%VARIABLES(14)
153 viscf=intbuf_tab%VARIABLES(15)
154 fnor =intbuf_tab%VARIABLES(4)
155 depth=intbuf_tab%VARIABLES(5)
156 fric_last =intbuf_tab%VARIABLES(6)
157 fnor_last =intbuf_tab%VARIABLES(7)
158C
159C
160 IF(nty==3)THEN
161C
162 ELSEIF(nty==8)THEN
163 ALLOCATE(
164 . distance(nsn),
165 . dist(nsn),
166 . ix1(nsn),
167 . ix2(nsn),
168 . ix3(nsn),
169 . ix4(nsn),
170 . x1(nsn), x2(nsn), x3(nsn), x4(nsn),
171 . y1(nsn), y2(nsn), y3(nsn), y4(nsn),
172 . z1(nsn), z2(nsn), z3(nsn), z4(nsn),
173 . xi(nsn), xi_a(nsn),
174 . yi(nsn), yi_a(nsn),
175 . zi(nsn), zi_a(nsn),
176 . xp(nsn), n1(nsn),
177 . yp(nsn), n2(nsn),
178 . zp(nsn), n3(nsn),
179 . ans(nsn),ssc(nsn),ttc(nsn),
180 . h1(nsn),h2(nsn),h3(nsn),h4(nsn),
181 . face(nsn),face_a(nsn),
182 . stif(nsn),
183 . x0(nsn), xx1(nsn), xx2(nsn),
184 . y0(nsn), yy1(nsn), yy2(nsn),
185 . z0(nsn), zz1(nsn), zz2(nsn),
186 . xx3(nsn), xx4(nsn), xi1(nsn),
187 . yy3(nsn), yy4(nsn), yi1(nsn),
188 . zz3(nsn), zz4(nsn), zi1(nsn),
189 . xi2(nsn), xi3(nsn), xi4(nsn),
190 . yi2(nsn), yi3(nsn), yi4(nsn),
191 . zi2(nsn), zi3(nsn), zi4(nsn),
192 . xn1(nsn), xn2(nsn), xn3(nsn), xn4(nsn),
193 . yn1(nsn), yn2(nsn), yn3(nsn), yn4(nsn),
194 . zn1(nsn), zn2(nsn), zn3(nsn), zn4(nsn),
195 . thk(nsn), area(nsn),alp(nsn),
196 . has_moved(nsn),
197 . tab_rmax(nsn),tab_rmax_uid(4,nsn),
198 . irtl_a(nsn),nsv_a(nsn),nsv2(nsn),
199 . index_in_comm(nmn),distlin(nsn))
200
201 dist = 0
202 distlin(1:nsn) = zero
203! DO 180 NG=1,NGROUS
204 nft=0
205 lft=1
206 llt=nsn
207 IF(nspmd > 1) THEN
208
209 CALL spmd_i8_index(nmn,t8%SPMD_COMM_PATTERN,index_in_comm,t8%S_COMM)
210
211 IF(t8%IS_ACTIVATED == 0) THEN
212 CALL spmd_i8_commslv(nsn,
213 . intbuf_tab%ILOCS,intbuf_tab%NSV,
214 . itab,t8%BUFFER,t8%SPMD_COMM_PATTERN,index_in_comm)
215
216 CALL spmd_i8_updbuf(nsn, intbuf_tab%ILOCS,intbuf_tab%NSV,
217 . itab,t8%BUFFER,t8%SPMD_COMM_PATTERN,index_in_comm)
218 t8%IS_ACTIVATED = 1
219 ENDIF
220
221 ENDIF
222
223 CALL i8loc3(
224 1 x, intbuf_tab%IRECTM,intbuf_tab%LMSR, intbuf_tab%MSR,
225 2 intbuf_tab%NSV, intbuf_tab%ILOCS, intbuf_tab%NSEGM, xi,
226 3 yi, zi, face, itab,
227 4 distance, iflinear, distlin, nsn,
228 5 lft, llt, nft)
229 IF(nspmd > 1) THEN
230
231 CALL spmd_i8_iloc(intbuf_tab%ILOCS, intbuf_tab%MSR,itab,t8%BUFFER,
232 . distance)
233
234 CALL spmd_i8_commslv(nsn,
235 . intbuf_tab%ILOCS,intbuf_tab%NSV,
236 . itab,t8%BUFFER,t8%SPMD_COMM_PATTERN,index_in_comm)
237
238 CALL spmd_i8_updbuf(nsn,intbuf_tab%ILOCS,intbuf_tab%NSV,
239 . itab,t8%BUFFER,t8%SPMD_COMM_PATTERN,index_in_comm)
240 ENDIF
241 CALL i8msr3(
242 1 x, intbuf_tab%IRECTM,intbuf_tab%LMSR, intbuf_tab%MSR,
243 2 intbuf_tab%NSV, intbuf_tab%ILOCS, intbuf_tab%IRTLM, intbuf_tab%NSEGM,
244 3 face, nsn, itab, has_moved,
245 4 tab_rmax, tab_rmax_uid, lft, llt,
246 5 nft)
247
248 IF(nspmd > 1) THEN
249 CALL spmd_i8_irtl(intbuf_tab%IRTLM,has_moved,
250 . tab_rmax,tab_rmax_uid,
251 . itab,t8%BUFFER)
252
253 ENDIF
254
255 nft=0
256 lft=1
257 llt=0
258 nsv_a(1:nsn) = 0
259 irtl_a(1:nsn) = 0
260 nsv2(1:nsn) = 0
261 xi_a(1:nsn) = 0
262 yi_a(1:nsn) = 0
263 zi_a(1:nsn) = 0
264 face_a(1:nsn) = 0
265
266 ! Compact NSV,IRTL,XI,ZI and FACE
267 ! by keeping only secnd nodes that are
268 ! active on this SPMD domain
269 DO i = 1,nsn
270 IF(intbuf_tab%ILOCS(i) > -1 .AND.
271 . intbuf_tab%IRTLM(i) > 0) THEN
272 llt = llt +1
273 nsv2(llt) = i
274 nsv_a(llt) = intbuf_tab%NSV(i)
275 irtl_a(llt) = intbuf_tab%IRTLM(i)
276 xi_a(llt) = xi(i)
277 yi_a(llt) = yi(i)
278 zi_a(llt) = zi(i)
279 face_a(llt) = face(i)
280
281 ELSEIF (iform==2) THEN
282 !WRITE(6,*) "S_FTSAVX=",INTBUF_TAB%S_FTSAVX,NSN
283 intbuf_tab%FTSAVX(i) = 0
284 intbuf_tab%FTSAVY(i) = 0
285 intbuf_tab%FTSAVZ(i) = 0
286 ENDIF
287 ENDDO
288
289 ninskid = 0
290 IF(h3d_data%N_SCAL_SKID > 0) THEN
291 ninskid = h3d_data%N_SKID_INTER(nin)
292 ENDIF
293
294 CALL i8cor3(
295 1 x, intbuf_tab%IRECTM,intbuf_tab%MSR, nsv_a,
296 2 irtl_a, ix1, ix2, ix3,
297 3 ix4, x1, y1, z1,
298 4 x2, y2, z2, x3,
299 5 y3, z3, x4, y4,
300 6 z4, lft, llt, nft)
301 CALL i8cst3(
302 1 x1, y1, z1, x2,
303 2 y2, z2, x3, y3,
304 3 z3, x4, y4, z4,
305 4 xi_a, yi_a, zi_a, n1,
306 5 n2, n3, ans, ssc,
307 6 ttc, face_a, x0, y0,
308 7 z0, xx1, yy1, zz1,
309 8 xx2, yy2, zz2, xx3,
310 9 yy3, zz3, xx4, yy4,
311 a zz4, xi1, yi1, zi1,
312 b xi2, yi2, zi2, xi3,
313 c yi3, zi3, xi4, yi4,
314 d zi4, xn1, yn1, zn1,
315 e xn2, yn2, zn2, xn3,
316 f yn3, zn3, xn4, yn4,
317 g zn4, area, lft, llt)
318
319 CALL i8gap3(
320 1 gap, thk, area, alp,
321 2 lft, llt)
322 CALL i8dis3(
323 1 igimp, nty, dist, x1,
324 2 y1, z1, x2, y2,
325 3 z2, x3, y3, z3,
326 4 x4, y4, z4, xi_a,
327 5 yi_a, zi_a, xp, yp,
328 6 zp, n1, n2, n3,
329 7 ans, ssc, ttc, h1,
330 8 h2, h3, h4, face_a,
331 9 alp, lft, llt)
332
333 CALL i8for3(lft ,llt ,nft ,a ,
334 2 intbuf_tab%MSR,nsv_a,irtl_a,intbuf_tab%STFM,
335 . intbuf_tab%NSV,nsv2, intbuf_tab%ILOCS,
336 3 intbuf_tab%STFNS,ibc ,icodt ,fsav ,igimp ,
337 4 x ,v ,ms ,fric ,nsn ,
338 5 fskyi ,isky ,fcont ,rcontact ,iform ,
339 6 intbuf_tab%FTSAVX,intbuf_tab%FTSAVY,intbuf_tab%FTSAVZ,visc ,fnor ,
340 7 depth ,dist ,intbuf_tab%GAPN,intbuf_tab%STF8,
341 8 stifn ,fncont ,ftcont ,itab, ift0,
342 9 ix1 ,ix2 ,ix3 ,ix4,
343 a xi_a ,yi_a ,zi_a,
344 b n1 ,n2 ,n3,
345 c ans , ssc , ttc,
346 d h1 , h2 , h3 , h4,
347 e face_a ,stif , xx3,
348 g yy3 , zz3 , xx4,
349 h yy4 , zz4 , xi1,
350 i yi1 , zi1 , xi2,
351 j yi2 , zi2 , xi3,
352 k yi3 , zi3 , xi4,
353 l thk ,h3d_data , ninskid,
354 m h3d_data%N_SCAL_SKID, pskids,intbuf_tab%IRECTM,nin,
355 n tagncont ,kloadpinter,loadpinter ,loadp_hyd_inter,
356 o iflinear ,fric_last,fnor_last,distlin)
357
358
359
360 DEALLOCATE(
361 . distance,
362 . dist,
363 . ix1,
364 . ix2,
365 . ix3,
366 . ix4,
367 . x1, x2, x3, x4,
368 . y1, y2, y3, y4,
369 . z1, z2, z3, z4,
370 . xi, xi_a,
371 . yi, yi_a,
372 . zi, zi_a,
373 . xp, n1,
374 . yp, n2,
375 . zp, n3,
376 . ans,ssc,ttc,
377 . h1,h2,h3,h4,
378 . face,face_a,
379 . stif,
380 . x0, xx1, xx2,
381 . y0, yy1, yy2,
382 . z0, zz1, zz2,
383 . xx3, xx4, xi1,
384 . yy3, yy4, yi1,
385 . zz3, zz4, zi1,
386 . xi2, xi3, xi4,
387 . yi2, yi3, yi4,
388 . zi2, zi3, zi4,
389 . xn1, xn2, xn3, xn4,
390 . yn1, yn2, yn3, yn4,
391 . zn1, zn2, zn3, zn4,
392 . thk, area,alp,
393 . has_moved,
394 . tab_rmax,tab_rmax_uid,
395 . irtl_a,nsv_a,nsv2,
396 . index_in_comm,
397 . distlin)
398
399 ENDIF ! TYPE 8
400
401 RETURN
402 END
#define my_real
Definition cppsort.cpp:32
subroutine area(d1, x, x2, y, y2, eint, stif0)
subroutine i8cor3(x, irect, msr, nsv, irtl, ix1, ix2, ix3, ix4, x1, y1, z1, x2, y2, z2, x3, y3, z3, x4, y4, z4, lft, llt, nft)
Definition i8cor3.F:35
subroutine i8cst3(x1, y1, z1, x2, y2, z2, x3, y3, z3, x4, y4, z4, xi, yi, zi, n1, n2, n3, ans, ssc, ttc, xface, x0, y0, z0, xx1, yy1, zz1, xx2, yy2, zz2, xx3, yy3, zz3, xx4, yy4, zz4, xi1, yi1, zi1, xi2, yi2, zi2, xi3, yi3, zi3, xi4, yi4, zi4, xn1, yn1, zn1, xn2, yn2, zn2, xn3, yn3, zn3, xn4, yn4, zn4, area, lft, llt)
Definition i8cst3.F:45
subroutine i8dis3(igimp, nty, dist, x1, y1, z1, x2, y2, z2, x3, y3, z3, x4, y4, z4, xi, yi, zi, xp, yp, zp, n1, n2, n3, ans, ssc, ttc, h1, h2, h3, h4, xface, alp, lft, llt)
Definition i8dis3.F:38
subroutine i8for3(lft, llt, nft, e, msr, nsv, irtl, stf, nsvglo, nsv2, iloc, stfn, ibc, icodt, fsav, igimp, x, v, ms, fmax, nsn, fskyi, isky, fcont, rcontact, iform, ftsavx, ftsavy, ftsavz, visc, fnor, depth, dist, gapn, slopen, stifn, fncont, ftcont, itab, ift0, ix1, ix2, ix3, ix4, xi, yi, zi, n1, n2, n3, ans, ssc, ttc, h1, h2, h3, h4, xface, stif, fni, fxi, fyi, fzi, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, thk, h3d_data, ninskid, ninterskid, pskids, irect, nin, tagncont, kloadpinter, loadpinter, loadp_hyd_inter, iflinear, fric_last, fnor_last, distlin)
Definition i8for3.F:57
subroutine i8gap3(gap, thk, area, alp, lft, llt)
Definition i8gap3.F:31
subroutine i8loc3(x, irect, lmsr, msr, nsv, iloc, nseg, xi, yi, zi, xface, itab, distance, iflinear, distlin, nsn, lft, llt, nft)
Definition i8loc3.F:36
subroutine i8msr3(x, irect, lmsr, msr, nsv, iloc, irtl, nseg, xface, nbsecnds, itab, has_moved, tab_rmax, tab_rmax_uid, lft, llt, nft)
Definition i8msr3.F:39
subroutine intvo8(ipari, x, a, icodt, fsav, v, ms, fskyi, isky, fcont, fncont, ftcont, icontact, rcontact, stifn, itab, intbuf_tab, t8, h3d_data, nin, pskids, tagncont, kloadpinter, loadpinter, loadp_hyd_inter)
Definition intvo8.F:52
subroutine spmd_i8_index(nmn, frontier, index_in_comm, s_comm)
subroutine spmd_i8_commslv(nbsecnds, iloc, nsv, itab, buffer, frontier, index_in_comm)
subroutine spmd_i8_updbuf(nbsecnds, iloc, nsv, itab, buffer, frontier, index_in_comm)
subroutine spmd_i8_iloc(iloc, msr, itab, buffer, distance)
Definition spmd_i8tool.F:34
subroutine spmd_i8_irtl(irtl, has_moved, tab_rmax, tab_rmax_uid, itab, buffer)