OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i10xsave.F File Reference
#include "implicit_f.inc"
#include "comlock.inc"
#include "com04_c.inc"
#include "task_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine i10xsave (x, nsv, msr, nsn, nmn, itask, xsav, xmin, ymin, zmin, xmax, ymax, zmax)

Function/Subroutine Documentation

◆ i10xsave()

subroutine i10xsave ( x,
integer, dimension(*) nsv,
integer, dimension(*) msr,
integer nsn,
integer nmn,
integer itask,
xsav,
xmin,
ymin,
zmin,
xmax,
ymax,
zmax )

Definition at line 30 of file i10xsave.F.

34C sauvegarde des XSAV et calcul borne domaine
35C-----------------------------------------------
36C I m p l i c i t T y p e s
37C-----------------------------------------------
38#include "implicit_f.inc"
39#include "comlock.inc"
40C-----------------------------------------------
41C C o m m o n B l o c k s
42C-----------------------------------------------
43#include "com04_c.inc"
44#include "task_c.inc"
45C-----------------------------------------------
46C D u m m y A r g u m e n t s
47C-----------------------------------------------
48 INTEGER NSN,NMN,ITASK,
49 . NSV(*),MSR(*)
51 . xmax, ymax, zmax, xmin, ymin, zmin,
52 . x(3,*), xsav(3,*)
53C-----------------------------------------------
54C L o c a l V a r i a b l e s
55C-----------------------------------------------
56 INTEGER NSNF,NMNF,NSNL,NMNL,I, J, II
58 . xxx, yyy, zzz
59C-----------------------------------------------
60C S o u r c e L i n e s
61C-----------------------------------------------
62C
63 nsnf=1+itask*nsn/nthread
64 nsnl=(itask+1)*nsn/nthread
65 nmnf=1+itask*nmn/nthread
66 nmnl=(itask+1)*nmn/nthread
67C
68 xmin=ep30
69 xmax=-ep30
70 ymin=ep30
71 ymax=-ep30
72 zmin=ep30
73 zmax=-ep30
74C
75 IF(nsn+nmn < numnod)THEN
76 DO i=nsnf,nsnl
77 j=nsv(i)
78 IF(j>0)THEN
79 xsav(1,i)=x(1,j)
80 xsav(2,i)=x(2,j)
81 xsav(3,i)=x(3,j)
82 END IF
83 END DO
84 DO i=nmnf,nmnl
85 ii = i+nsn
86 j=msr(i)
87C msr < 0 <=> shooting nodes
88 IF(j>0) THEN
89 xmin= min(xmin,x(1,j))
90 ymin= min(ymin,x(2,j))
91 zmin= min(zmin,x(3,j))
92 xmax= max(xmax,x(1,j))
93 ymax= max(ymax,x(2,j))
94 zmax= max(zmax,x(3,j))
95 xsav(1,ii)=x(1,j)
96 xsav(2,ii)=x(2,j)
97 xsav(3,ii)=x(3,j)
98 ENDIF
99 END DO
100 ELSE
101C
102 DO i=nsnf,nsnl
103 j=nsv(i)
104 IF(j>0)THEN
105 xsav(1,j)=x(1,j)
106 xsav(2,j)=x(2,j)
107 xsav(3,j)=x(3,j)
108 END IF
109 END DO
110 DO i=nmnf,nmnl
111 j=msr(i)
112C msr < 0 <=> shooting nodes
113 IF(j>0) THEN
114 xmin= min(xmin,x(1,j))
115 ymin= min(ymin,x(2,j))
116 zmin= min(zmin,x(3,j))
117 xmax= max(xmax,x(1,j))
118 ymax= max(ymax,x(2,j))
119 zmax= max(zmax,x(3,j))
120 xsav(1,j)=x(1,j)
121 xsav(2,j)=x(2,j)
122 xsav(3,j)=x(3,j)
123 END IF
124 END DO
125 ENDIF
126C
127 RETURN
#define my_real
Definition cppsort.cpp:32
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