OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
s20mass3.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "vect01_c.inc"
#include "param_c.inc"
#include "com01_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine s20mass3 (mass, ms, partsav, ipart, mss, volg, xx, yy, zz, vx, vy, vz, nc, sti, stifn, deltax2, rho, dtx, dtelem, mssx, rhocp, mcp, mcps, mcpsx, fill)
subroutine s20msi (rho, mass, volu, dtelem, sti, off, sig, eint, dtx, nel, offg, sigg, eintg, rhog, wip)
subroutine sigin20b (sig, pm, vol, sigsp, sigi, eint, rho, uvar, eps, ix, nix, nsigi, ipt, nuvar, nel, iuser, idef, nsigs, strsglob, straglob, jhbe, igtyp, x, bufgama, mat, epsp, l_pla, pt, sigb, l_sigb, ipm, bufmat, voldp)

Function/Subroutine Documentation

◆ s20mass3()

subroutine s20mass3 ( mass,
ms,
partsav,
integer, dimension(*) ipart,
mss,
volg,
xx,
yy,
zz,
vx,
vy,
vz,
integer, dimension(mvsiz,20) nc,
sti,
stifn,
deltax2,
rho,
dtx,
dtelem,
mssx,
rhocp,
mcp,
mcps,
mcpsx,
fill )

Definition at line 28 of file s20mass3.F.

34C-----------------------------------------------
35C I m p l i c i t T y p e s
36C-----------------------------------------------
37#include "implicit_f.inc"
38C-----------------------------------------------
39C G l o b a l P a r a m e t e r s
40C-----------------------------------------------
41#include "mvsiz_p.inc"
42C-----------------------------------------------
43C D u m m y A r g u m e n t s
44C-----------------------------------------------
45 INTEGER IPART(*), NC(MVSIZ,20)
47 . mass(*), ms(*),partsav(20,*), mss(8,*),deltax2(*),
48 . xx(mvsiz,20), yy(mvsiz,20), zz(mvsiz,20),
49 . vx(mvsiz,20), vy(mvsiz,20), vz(mvsiz,20),sti(*),stifn(*),
50 . volg(mvsiz),rho(mvsiz),dtx(mvsiz),dtelem(mvsiz),mssx(12,*),
51 . rhocp(*), mcp(*), mcps(8,*),mcpsx(12,*), fill(*)
52C-----------------------------------------------
53C C o m m o n B l o c k s
54C-----------------------------------------------
55#include "vect01_c.inc"
56C-----------------------------------------------
57C L o c a l V a r i a b l e s
58C-----------------------------------------------
59 INTEGER I, IP,N1,N2, IPERM1(20),IPERM2(20),N
60C REAL
62 . axx,ayy,azz,axy,ayz,azx,am,bm,fac,masscp
63C
64 DATA iperm1/0,0,0,0,0,0,0,0,1,2,3,4,1,2,3,4,5,6,7,8/
65 DATA iperm2/0,0,0,0,0,0,0,0,2,3,4,1,5,6,7,8,6,7,8,5/
66C-----------------------------------------------------------------------
67C mass init in parith for spmd
68 DO i=lft,llt
69C
70 mass(i)=fill(i)*rho(i)*volg(i)
71 dtelem(i)=min(dtelem(i),dtx(i))
72 sti(i) = fill(i) * rho(i) * volg(i) * two / sixty4 /
73 . max(em20,dtx(i)*dtx(i))
74 fac=three_over_14
75 am = mass(i)*fac/(eight*fac + twelve)
76 bm = mass(i)*one/(eight*fac + twelve)
77 mss(1,i)=am
78 mss(2,i)=am
79 mss(3,i)=am
80 mss(4,i)=am
81 mss(5,i)=am
82 mss(6,i)=am
83 mss(7,i)=am
84 mss(8,i)=am
85 stifn(nc(i,1))=stifn(nc(i,1))+sti(i)*deltax2(i)
86 stifn(nc(i,2))=stifn(nc(i,2))+sti(i)*deltax2(i)
87 stifn(nc(i,3))=stifn(nc(i,3))+sti(i)*deltax2(i)
88 stifn(nc(i,4))=stifn(nc(i,4))+sti(i)*deltax2(i)
89 stifn(nc(i,5))=stifn(nc(i,5))+sti(i)*deltax2(i)
90 stifn(nc(i,6))=stifn(nc(i,6))+sti(i)*deltax2(i)
91 stifn(nc(i,7))=stifn(nc(i,7))+sti(i)*deltax2(i)
92 stifn(nc(i,8))=stifn(nc(i,8))+sti(i)*deltax2(i)
93 DO n=9,20
94 n1=iperm1(n)
95 n2=iperm2(n)
96 IF(nc(i,n) /= 0)THEN
97 mssx(n-8,i)=bm
98 stifn(nc(i,n))=stifn(nc(i,n))+sti(i)
99 ELSE
100 mss(n1,i)=mss(n1,i)+ half*bm
101 mss(n2,i)=mss(n2,i)+ half*bm
102 stifn(nc(i,n1))=stifn(nc(i,n1)) + half*sti(i)
103 stifn(nc(i,n2))=stifn(nc(i,n2)) + half*sti(i)
104 ENDIF
105 ENDDO
106C
107 ip=ipart(i)
108 partsav(1,ip)=partsav(1,ip) + mass(i)
109 partsav(2,ip)=partsav(2,ip)
110 . + am*(xx(i,1)+xx(i,2)+xx(i,3)+xx(i,4)
111 . +xx(i,5)+xx(i,6)+xx(i,7)+xx(i,8))
112 . + bm*(xx(i,9) +xx(i,10)+xx(i,11)+xx(i,12)+xx(i,13)+xx(i,14)
113 . +xx(i,15)+xx(i,16)+xx(i,17)+xx(i,18)+xx(i,19)+xx(i,20))
114 partsav(3,ip)=partsav(3,ip)
115 . + am*(yy(i,1)+yy(i,2)+yy(i,3)+yy(i,4)
116 . +yy(i,5)+yy(i,6)+yy(i,7)+yy(i,8))
117 . + bm*(yy(i,9) +yy(i,10)+yy(i,11)+yy(i,12)+yy(i,13)+yy(i,14)
118 . +yy(i,15)+yy(i,16)+yy(i,17)+yy(i,18)+yy(i,19)+yy(i,20))
119 partsav(4,ip)=partsav(4,ip)
120 . + am*(zz(i,1)+zz(i,2)+zz(i,3)+zz(i,4)
121 . +zz(i,5)+zz(i,6)+zz(i,7)+zz(i,8))
122 . + bm*(zz(i,9) +zz(i,10)+zz(i,11)+zz(i,12)+zz(i,13)+zz(i,14)
123 . +zz(i,15)+zz(i,16)+zz(i,17)+zz(i,18)+zz(i,19)+zz(i,20))
124 axx = am*(xx(i,1)*xx(i,1)+xx(i,2)*xx(i,2)
125 . +xx(i,3)*xx(i,3)+xx(i,4)*xx(i,4)
126 . +xx(i,5)*xx(i,5)+xx(i,6)*xx(i,6)
127 . +xx(i,7)*xx(i,7)+xx(i,8)*xx(i,8))
128 . +bm*(xx(i,9) *xx(i,9) +xx(i,10)*xx(i,10)
129 . +xx(i,11)*xx(i,11)+xx(i,12)*xx(i,12)
130 . +xx(i,13)*xx(i,13)+xx(i,14)*xx(i,14)
131 . +xx(i,15)*xx(i,15)+xx(i,16)*xx(i,16)
132 . +xx(i,17)*xx(i,17)+xx(i,18)*xx(i,18)
133 . +xx(i,19)*xx(i,19)+xx(i,20)*xx(i,20))
134 ayy = am*(yy(i,1)*yy(i,1)+yy(i,2)*yy(i,2)
135 . +yy(i,3)*yy(i,3)+yy(i,4)*yy(i,4)
136 . +yy(i,5)*yy(i,5)+yy(i,6)*yy(i,6)
137 . +yy(i,7)*yy(i,7)+yy(i,8)*yy(i,8))
138 . +bm*(yy(i,9) *yy(i,9) +yy(i,10)*yy(i,10)
139 . +yy(i,11)*yy(i,11)+yy(i,12)*yy(i,12)
140 . +yy(i,13)*yy(i,13)+yy(i,14)*yy(i,14)
141 . +yy(i,15)*yy(i,15)+yy(i,16)*yy(i,16)
142 . +yy(i,17)*yy(i,17)+yy(i,18)*yy(i,18)
143 . +yy(i,19)*yy(i,19)+yy(i,20)*yy(i,20))
144 azz = am*(zz(i,1)*zz(i,1)+zz(i,2)*zz(i,2)
145 . +zz(i,3)*zz(i,3)+zz(i,4)*zz(i,4)
146 . +zz(i,5)*zz(i,5)+zz(i,6)*zz(i,6)
147 . +zz(i,7)*zz(i,7)+zz(i,8)*zz(i,8))
148 . +bm*(zz(i,9) *zz(i,9) +zz(i,10)*zz(i,10)
149 . +zz(i,11)*zz(i,11)+zz(i,12)*zz(i,12)
150 . +zz(i,13)*zz(i,13)+zz(i,14)*zz(i,14)
151 . +zz(i,15)*zz(i,15)+zz(i,16)*zz(i,16)
152 . +zz(i,17)*zz(i,17)+zz(i,18)*zz(i,18)
153 . +zz(i,19)*zz(i,19)+zz(i,20)*zz(i,20))
154 axy = am*(xx(i,1)*yy(i,1)+xx(i,2)*yy(i,2)
155 . +xx(i,3)*yy(i,3)+xx(i,4)*yy(i,4)
156 . +xx(i,5)*yy(i,5)+xx(i,6)*yy(i,6)
157 . +xx(i,7)*yy(i,7)+xx(i,8)*yy(i,8))
158 . +bm*(xx(i,9) *yy(i,9) +xx(i,10)*yy(i,10)
159 . +xx(i,11)*yy(i,11)+xx(i,12)*yy(i,12)
160 . +xx(i,13)*yy(i,13)+xx(i,14)*yy(i,14)
161 . +xx(i,15)*yy(i,15)+xx(i,16)*yy(i,16)
162 . +xx(i,17)*yy(i,17)+xx(i,18)*yy(i,18)
163 . +xx(i,19)*yy(i,19)+xx(i,20)*yy(i,20))
164 ayz = am*(yy(i,1)*zz(i,1)+yy(i,2)*zz(i,2)
165 . +yy(i,3)*zz(i,3)+yy(i,4)*zz(i,4)
166 . +yy(i,5)*zz(i,5)+yy(i,6)*zz(i,6)
167 . +yy(i,7)*zz(i,7)+yy(i,8)*zz(i,8))
168 . +bm*(yy(i,9) *zz(i,9) +yy(i,10)*zz(i,10)
169 . +yy(i,11)*zz(i,11)+yy(i,12)*zz(i,12)
170 . +yy(i,13)*zz(i,13)+yy(i,14)*zz(i,14)
171 . +yy(i,15)*zz(i,15)+yy(i,16)*zz(i,16)
172 . +yy(i,17)*zz(i,17)+yy(i,18)*zz(i,18)
173 . +yy(i,19)*zz(i,19)+yy(i,20)*zz(i,20))
174 azx = am*(zz(i,1)*xx(i,1)+zz(i,2)*xx(i,2)
175 . +zz(i,3)*xx(i,3)+zz(i,4)*xx(i,4)
176 . +zz(i,5)*xx(i,5)+zz(i,6)*xx(i,6)
177 . +zz(i,7)*xx(i,7)+zz(i,8)*xx(i,8))
178 . +bm*(zz(i,9) *xx(i,9) +zz(i,10)*xx(i,10)
179 . +zz(i,11)*xx(i,11)+zz(i,12)*xx(i,12)
180 . +zz(i,13)*xx(i,13)+zz(i,14)*xx(i,14)
181 . +zz(i,15)*xx(i,15)+zz(i,16)*xx(i,16)
182 . +zz(i,17)*xx(i,17)+zz(i,18)*xx(i,18)
183 . +zz(i,19)*xx(i,19)+zz(i,20)*xx(i,20))
184 partsav(5,ip) =partsav(5,ip) + (ayy+azz)
185 partsav(6,ip) =partsav(6,ip) + (azz+axx)
186 partsav(7,ip) =partsav(7,ip) + (axx+ayy)
187 partsav(8,ip) =partsav(8,ip) - axy
188 partsav(9,ip) =partsav(9,ip) - ayz
189 partsav(10,ip)=partsav(10,ip) - azx
190C
191 partsav(11,ip)=partsav(11,ip)
192 . + am*(vx(i,1)+vx(i,2)+vx(i,3)+vx(i,4)
193 . +vx(i,5)+vx(i,6)+vx(i,7)+vx(i,8))
194 . + bm*(vx(i,9) +vx(i,10)+vx(i,11)+vx(i,12)+vx(i,13)+vx(i,14)
195 . +vx(i,15)+vx(i,16)+vx(i,17)+vx(i,18)+vx(i,19)+vx(i,20))
196 partsav(12,ip)=partsav(12,ip)
197 . + am*(vy(i,1)+vy(i,2)+vy(i,3)+vy(i,4)
198 . +vy(i,5)+vy(i,6)+vy(i,7)+vy(i,8))
199 . + bm*(vy(i,9) +vy(i,10)+vy(i,11)+vy(i,12)+vy(i,13)+vy(i,14)
200 . +vy(i,15)+vy(i,16)+vy(i,17)+vy(i,18)+vy(i,19)+vy(i,20))
201 partsav(13,ip)=partsav(13,ip)
202 . + am*(vz(i,1)+vz(i,2)+vz(i,3)+vz(i,4)
203 . +vz(i,5)+vz(i,6)+vz(i,7)+vz(i,8))
204 . + bm*(vz(i,9) +vz(i,10)+vz(i,11)+vz(i,12)+vz(i,13)+vz(i,14)
205 . +vz(i,15)+vz(i,16)+vz(i,17)+vz(i,18)+vz(i,19)+vz(i,20))
206 partsav(14,ip)=partsav(14,ip) + half *
207 . (am*(vx(i,1)*vx(i,1)+vx(i,2)*vx(i,2)
208 . +vx(i,3)*vx(i,3)+vx(i,4)*vx(i,4)
209 . +vx(i,5)*vx(i,5)+vx(i,6)*vx(i,6)
210 . +vx(i,7)*vx(i,7)+vx(i,8)*vx(i,8)
211 . +vy(i,1)*vy(i,1)+vy(i,2)*vy(i,2)
212 . +vy(i,3)*vy(i,3)+vy(i,4)*vy(i,4)
213 . +vy(i,5)*vy(i,5)+vy(i,6)*vy(i,6)
214 . +vy(i,7)*vy(i,7)+vy(i,8)*vy(i,8)
215 . +vz(i,1)*vz(i,1)+vz(i,2)*vz(i,2)
216 . +vz(i,3)*vz(i,3)+vz(i,4)*vz(i,4)
217 . +vz(i,5)*vz(i,5)+vz(i,6)*vz(i,6)
218 . +vz(i,7)*vz(i,7)+vz(i,8)*vz(i,8))
219 . +bm*(vx(i,9) *vx(i,9) +vx(i,10)*vx(i,10)
220 . +vx(i,11)*vx(i,11)+vx(i,12)*vx(i,12)
221 . +vx(i,13)*vx(i,13)+vx(i,14)*vx(i,14)
222 . +vx(i,15)*vx(i,15)+vx(i,16)*vx(i,16)
223 . +vx(i,17)*vx(i,17)+vx(i,18)*vx(i,18)
224 . +vx(i,19)*vx(i,19)+vx(i,20)*vx(i,20)
225 . +vy(i,9) *vy(i,9) +vy(i,10)*vy(i,10)
226 . +vy(i,11)*vy(i,11)+vy(i,12)*vy(i,12)
227 . +vy(i,13)*vy(i,13)+vy(i,14)*vy(i,14)
228 . +vy(i,15)*vy(i,15)+vy(i,16)*vy(i,16)
229 . +vy(i,17)*vy(i,17)+vy(i,18)*vy(i,18)
230 . +vy(i,19)*vy(i,19)+vy(i,20)*vy(i,20)
231 . +vz(i,9) *vz(i,9) +vz(i,10)*vz(i,10)
232 . +vz(i,11)*vz(i,11)+vz(i,12)*vz(i,12)
233 . +vz(i,13)*vz(i,13)+vz(i,14)*vz(i,14)
234 . +vz(i,15)*vz(i,15)+vz(i,16)*vz(i,16)
235 . +vz(i,17)*vz(i,17)+vz(i,18)*vz(i,18)
236 . +vz(i,19)*vz(i,19)+vz(i,20)*vz(i,20)))
237 ENDDO
238C
239c --- for FEM solide heat trasnfert
240c
241 IF(jthe < 0 ) THEN
242 DO i=lft,llt
243 masscp=fill(i)*rhocp(i)*volg(i)
244 fac=three_over_14
245 am = masscp*fac/(eight*fac + twelve)
246 bm = masscp*one/(eight*fac + twelve)
247 mcps(1,i)=am
248 mcps(2,i)=am
249 mcps(3,i)=am
250 mcps(4,i)=am
251 mcps(5,i)=am
252 mcps(6,i)=am
253 mcps(7,i)=am
254 mcps(8,i)=am
255 DO n=9,20
256 n1=iperm1(n)
257 n2=iperm2(n)
258 IF(nc(i,n) /= 0)THEN
259 mcpsx(n-8,i)=bm
260 ELSE
261 mcps(n1,i)=mcps(n1,i)+ half*bm
262 mcps(n2,i)=mcps(n2,i)+ half*bm
263 ENDIF
264 ENDDO
265 ENDDO
266 ENDIF
267C-----------
268 RETURN
#define my_real
Definition cppsort.cpp:32
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21

◆ s20msi()

subroutine s20msi ( rho,
mass,
volu,
dtelem,
sti,
off,
sig,
eint,
dtx,
integer nel,
offg,
sigg,
eintg,
rhog,
wip )

Definition at line 278 of file s20mass3.F.

281C
282C-----------------------------------------------
283C I m p l i c i t T y p e s
284C-----------------------------------------------
285#include "implicit_f.inc"
286C-----------------------------------------------
287C G l o b a l P a r a m e t e r s
288C-----------------------------------------------
289#include "mvsiz_p.inc"
290C-----------------------------------------------
291C D u m m y A r g u m e n t s
292C-----------------------------------------------
293 INTEGER NEL
294C REAL
295 my_real
296 . rho(*), mass(*),volu(*),dtelem(*),sti(*),
297 . sig(nel,6),eint(*),off(*),
298 . sigg(nel,6),eintg(*),rhog(*),offg(*),wip
299C-----------------------------------------------
300C C o m m o n B l o c k s
301C-----------------------------------------------
302#include "vect01_c.inc"
303C-----------------------------------------------
304C L o c a l V a r i a b l e s
305C-----------------------------------------------
306 INTEGER I
307 my_real
308 . dtx(mvsiz)
309C-----------------------------------------------
310 DO i=lft,llt
311c MASS(I)=RHO(I)*VOLU(I)
312c DTELEM(I)=MIN(DTELEM(I),DTX(I))
313c STI(I) = RHO(I) * VOLU(I) /
314c . MAX(EM20,DTX(I)*DTX(I))
315 sigg(i,1) = sigg(i,1) + wip * sig(i,1)
316 sigg(i,2) = sigg(i,2) + wip * sig(i,2)
317 sigg(i,3) = sigg(i,3) + wip * sig(i,3)
318 sigg(i,4) = sigg(i,4) + wip * sig(i,4)
319 sigg(i,5) = sigg(i,5) + wip * sig(i,5)
320 sigg(i,6) = sigg(i,6) + wip * sig(i,6)
321 rhog(i) = rhog(i) + wip * rho(i)
322 eintg(i) = eintg(i) + wip * eint(i)
323 offg(i) = off(i)
324 ENDDO
325C-----------
326 RETURN

◆ sigin20b()

subroutine sigin20b ( sig,
pm,
vol,
sigsp,
sigi,
eint,
rho,
uvar,
eps,
integer, dimension(nix,*) ix,
integer nix,
integer nsigi,
integer ipt,
integer nuvar,
integer nel,
integer iuser,
integer idef,
integer nsigs,
integer, dimension(*) strsglob,
integer, dimension(*) straglob,
integer jhbe,
integer igtyp,
x,
bufgama,
integer, dimension(nel) mat,
epsp,
integer l_pla,
integer, dimension(*) pt,
integer, dimension(nel*l_sigb) sigb,
integer l_sigb,
integer, dimension(npropmi,*) ipm,
bufmat,
double precision, dimension(*) voldp )

Definition at line 343 of file s20mass3.F.

351C-----------------------------------------------
352C I m p l i c i t T y p e s
353C-----------------------------------------------
354#include "implicit_f.inc"
355C-----------------------------------------------
356C C o m m o n B l o c k s
357C-----------------------------------------------
358#include "param_c.inc"
359#include "com01_c.inc"
360#include "vect01_c.inc"
361C-----------------------------------------------
362C D u m m y A r g u m e n t s
363C-----------------------------------------------
364 INTEGER IPM(NPROPMI,*)
365 INTEGER NIX,JPS,JHBE,IGTYP,NSIGI,NUVAR,NEL,IUSER,NSIGS,
366 . IDEF,L_PLA,L_SIGB
367 INTEGER IX(NIX,*),STRSGLOB(*),STRAGLOB(*),MAT(NEL),PT(*),SIGB(NEL*L_SIGB)
368C REAL
369 my_real
370 . sig(nel,6) , eint(*), rho(*), sigsp(nsigi,*),
371 . pm(npropm,*), vol(*), uvar(*), sigi(nsigs,*),
372 . eps(nel,6),x(3,*),bufgama(*),epsp(*),bufmat(*)
373 double precision
374 . voldp(*)
375C-----------------------------------------------
376C L o c a l V a r i a b l e s
377C-----------------------------------------------
378 INTEGER I, II, JJ, IPT, IPP,IUS,IPSU,JPS1,MA,IFLAGINI,
379 . NVAR_TMP,IADB,NRATE
380 my_real
381 . gama(6),tens(6)
382C=======================================================================
383 DO i=lft,llt
384 ma = mat(i)
385 eint(i)=pm(23,ma)
386 rho(i) =pm(89,ma)
387 ENDDO
388C
389 IF (isigi /= 0) THEN
390C
391 jps=1+ (ipt-1)*9
392 jps1 = nvsolid1 + (ipt-1)*6
393 DO i=lft,llt
394 iflagini = 0
395 IF (straglob(i) == 1 .OR. strsglob(i) == 1) THEN
396 IF(jcvt==2)THEN
397 gama(1)=bufgama(i )
398 gama(2)=bufgama(i + nel)
399 gama(3)=bufgama(i + 2*nel)
400 gama(4)=bufgama(i + 3*nel)
401 gama(5)=bufgama(i + 4*nel)
402 gama(6)=bufgama(i + 5*nel)
403 ELSE
404 gama(1)=one
405 gama(2)=zero
406 gama(3)=zero
407 gama(4)=zero
408 gama(5)=one
409 gama(6)=zero
410 END IF
411 ENDIF
412C CONTRAINTES INITIALES
413
414 ii=nft+i
415 jj=pt(ii)
416 iflagini = 1
417 IF(jj==0)iflagini = 0
418c
419 IF (iflagini == 1) THEN
420 ipp = i
421 IF(nvsolid1 /= 0 )THEN
422 sig(i,1)=sigsp(jps+1,jj)
423 sig(i,2)=sigsp(jps+2,jj)
424 sig(i,3)=sigsp(jps+3,jj)
425 sig(i,4)=sigsp(jps+4,jj)
426 sig(i,5)=sigsp(jps+5,jj)
427 sig(i,6)=sigsp(jps+6,jj)
428 IF (strsglob(i) == 1) THEN
429 tens(1)=sig(i,1)
430 tens(2)=sig(i,2)
431 tens(3)=sig(i,3)
432 tens(4)=sig(i,4)
433 tens(5)=sig(i,5)
434 tens(6)=sig(i,6)
435 CALL srota6_m1(x,ix(1,ii),jcvt,tens,gama,jhbe,igtyp)
436 sig(i,1)=tens(1)
437 sig(i,2)=tens(2)
438 sig(i,3)=tens(3)
439 sig(i,4)=tens(4)
440 sig(i,5)=tens(5)
441 sig(i,6)=tens(6)
442 ENDIF
443 IF(l_pla /= 0 .AND. sigsp(jps+7,jj) /= zero)
444 . epsp(i) = sigsp(jps+7,jj)
445 IF (sigsp(jps+8,jj) /= zero) eint(i)=sigsp(jps+8,jj)
446 IF (sigsp(jps+9,jj) /= zero) THEN
447 voldp(i) = sigsp(jps+9,jj)*voldp(i) / rho(i)
448 vol(i) = sigsp(jps+9,jj)*vol(i) / rho(i)
449 rho(i) = sigsp(jps+9,jj)
450 ENDIF
451 ENDIF
452c
453 IF (mtn >= 28 .AND. iuser == 1) THEN
454 IF (mtn == 36 .and. l_sigb == 6) THEN
455 iadb = ipm(7,mat(1))
456 nrate = nint(bufmat(iadb))
457 nvar_tmp = sigsp(nvsolid1 + nvsolid2 + 3, jj)
458 ipsu = nvsolid1 + nvsolid2 + 4 + (ipt - 1)*nvar_tmp
459 DO ius = 1,6
460 ipp = i + (ius -1)*nel
461 sigb(ipp) = sigsp(ipsu + nrate + ius, jj)
462 ENDDO
463 ELSEIF (mtn == 112) THEN
464 nvar_tmp = sigsp(nvsolid1 + nvsolid2 + 3, jj)
465 ipsu = nvsolid1 + nvsolid2 + 4 + (ipt - 1)*nvar_tmp
466 DO ius = 1, nvar_tmp
467 ipp = i + ius*nel
468 epsp(ipp) = sigsp(ipsu + ius, jj)
469 ENDDO
470 ELSE IF (mtn /= 36) THEN
471 nvar_tmp = sigsp(nvsolid1 + nvsolid2 + 3, jj)
472 ipsu = nvsolid1 + nvsolid2 + 4 + (ipt - 1)*nvar_tmp
473 DO ius = 1, nvar_tmp
474 ipp = i + (ius -1)*nel
475 uvar(ipp) = sigsp(ipsu + ius, jj)
476 ENDDO
477 DO ius = nvar_tmp + 1, nuvar
478 ipp = i + (ius -1)*nel
479 uvar(ipp) = zero
480 ENDDO
481 ENDIF
482 ENDIF
483c
484 IF (nvsolid2 /= 0 .AND. idef /=0) THEN
485 eps(i,1) = sigsp(jps1 + 1 , jj)
486 eps(i,2) = sigsp(jps1 + 2 , jj)
487 eps(i,3) = sigsp(jps1 + 3 , jj)
488 eps(i,4) = sigsp(jps1 + 4 , jj)
489 eps(i,5) = sigsp(jps1 + 5 , jj)
490 eps(i,6) = sigsp(jps1 + 6 , jj)
491 IF (straglob(i) == 1) THEN
492 tens(1)=eps(i,1)
493 tens(2)=eps(i,2)
494 tens(3)=eps(i,3)
495 tens(4)=eps(i,4)
496 tens(5)=eps(i,5)
497 tens(6)=eps(i,6)
498 CALL srota6_m1(x,ix(1,ii),jcvt,tens,gama,jhbe,igtyp)
499 eps(i,1) = tens(1)
500 eps(i,2) = tens(2)
501 eps(i,3) = tens(3)
502 eps(i,4) = tens(4)
503 eps(i,5) = tens(5)
504 eps(i,6) = tens(6)
505 ENDIF
506 ENDIF
507 ENDIF !IF (IFLAGINI == 1)
508c
509 ENDDO ! I=LFT,LLT
510 ENDIF ! ISIGI /= 0
511c-----------
512 RETURN
subroutine srota6_m1(x, ixs, kcvt, tens, gama, khbe, ityp)
Definition srota6_M1.F:37