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