OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
volpvg.F File Reference
#include "implicit_f.inc"
#include "param_c.inc"
#include "com06_c.inc"
#include "com08_c.inc"
#include "task_c.inc"
#include "units_c.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "vect01_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine volpvga (ivolu, rvolu, vol, fsav, nvent, ibaghol, rbaghol, pmain, wfext)
subroutine volpvgb (ivolu, rvolu, vol, fsav, nvent, ibaghol, rbaghol, normal, nn, igrsurf, iparg, elbuf_tab, fr_mv, igroupc, igrouptg)

Function/Subroutine Documentation

◆ volpvga()

subroutine volpvga ( integer, dimension(*) ivolu,
rvolu,
vol,
fsav,
integer nvent,
integer, dimension(nibhol,*) ibaghol,
rbaghol,
integer pmain,
double precision, intent(inout) wfext )

Definition at line 29 of file volpvg.F.

31C-----------------------------------------------
32C MONITORED VOLUMES TYPE PERFECT GAS, INPUT STARTER >= 4.4 (PRESSURE).
33C-----------------------------------------------
34C I m p l i c i t T y p e s
35C-----------------------------------------------
36#include "implicit_f.inc"
37C-----------------------------------------------
38C C o m m o n B l o c k s
39C-----------------------------------------------
40#include "param_c.inc"
41#include "com06_c.inc"
42#include "com08_c.inc"
43#include "task_c.inc"
44C-----------------------------------------------
45C D u m m y A r g u m e n t s
46C-----------------------------------------------
47 INTEGER IVOLU(*),NVENT,IBAGHOL(NIBHOL,*),PMAIN
48 my_real rvolu(*),vol,fsav(*),rbaghol(nrbhol,*)
49C-----------------------------------------------
50C L o c a l V a r i a b l e s
51C-----------------------------------------------
52 INTEGER IDEF,IV, II
53 my_real vinc,gama,pres,pmax,veps,area,pold,vold,
54 . q,qold,dv,amu,rot,
55 . amtot,energy,energ_old,dmout,deout,fac,
56 . pdef,dtpdefc,m0, vol0
57 my_real trelax,dein, pext, pini, cv, temperature
58 INTEGER :: IEQUI
59 my_real :: min_tvent
60 DOUBLE PRECISION,INTENT(INOUT) :: WFEXT
61C-----------------------------------------------
62 pres = -huge(pres)
63 min_tvent = ep30
64 DO ii = 1, nvent
65 min_tvent = min(min_tvent, rbaghol(3, ii))
66 ENDDO
67
68 pext = rvolu(3)
69 pini = rvolu(56)
70 pold =rvolu(12)
71 m0 = rvolu(57)
72 vol0 = rvolu(4)
73 DO iv=1,nvent
74 idef = ibaghol(1,iv)
75 pdef = rbaghol(1,iv)
76 dtpdefc= rbaghol(5,iv)
77 IF(idef <= 0 .AND. pold > pdef + pext)
78 . rbaghol(5,iv)=dtpdefc+dt1
79 ENDDO
80C------------------------
81 idef =ivolu(14)
82 iequi = ivolu(15)
83 gama =rvolu(1)
84 vinc =rvolu(5)
85 pmax =rvolu(6)
86 energ_old=rvolu(13)
87 vold =rvolu(16)
88 veps =rvolu(17)
89 area =rvolu(18)
90 qold =rvolu(23)
91 vol =vol + veps
92 dmout =rvolu(24)
93 deout =rvolu(22)
94 amtot =rvolu(20)
95 amu =rvolu(2)
96 rot =rvolu(21)
97 trelax =rvolu(48)
98 cv = rvolu(19)
99 dv =vol-vold
100 temperature = zero
101 energy = zero
102C
103 IF(idef==1)THEN
104 pres = pext
105 q = zero
106 ELSE
107 IF (iequi == 0 .OR. trelax == zero) THEN
108 amtot=amtot-dmout*dt1
109C
110C CALCUL DE L ENERGIE PUIS DE LA PRESSION
111 fac = half*(gama-one)*dv
112 IF(trelax == zero .OR. tt > trelax)THEN
113 dein = zero
114 ELSE
115 dein = rvolu(52)
116 ENDIF
117 energy= ((one-fac/(vold-vinc))*energ_old +(dein-deout)*dt1) /
118 . (one+fac/(vol-vinc))
119 energy = max(energy,zero)
120C
121 pres=(gama-one)*energy/(vol-vinc)
122 ELSE IF (iequi == 1) THEN
123 IF (tt <= trelax) THEN
124 pres = pext + tt * (pini - pext) / trelax
125 energy = pres * (vol - vinc) / (gama - one)
126 amtot = m0 * pres * vol / (pext * vol0)
127 ELSE
128 pres = pext * amtot / m0 * vol0 / vol
129 energy = pres * (vol - vinc) / (gama - one)
130 ENDIF
131 ELSE IF (iequi == 2) THEN
132 IF (tt <= trelax) THEN
133 pres = pext + tt * (pini - pext) / trelax
134 energy = pres * (vol - vinc) / (gama - one)
135 amtot = m0 * (pres / pext)**(one / gama) * vol / vol0
136 ELSE
137 pres = pext * (amtot / m0 * vol0 / vol)**gama
138 energy = pres * (vol - vinc) / (gama - one)
139 ENDIF
140 ENDIF
141 IF (iequi /= 0) THEN
142 IF (tt > min_tvent .AND. tt > trelax) THEN
143! IEQUI = 0
144 ivolu(15) = 0
145 ENDIF
146 ENDIF
147
148 IF(dt1==zero.OR.dv>zero)THEN
149 q=zero
150 ELSE
151 q=-amu*sqrt(pres*area*rot/vol)*dv/area/dt1
152 ENDIF
153C
154 IF(pres>pmax)THEN
155 idef=1
156 pres = pext
157 q = zero
158 ENDIF
159 IF (iequi > 0) THEN
160 temperature = energy / (cv * amtot)
161 ENDIF
162 rvolu(20) = amtot
163 rvolu(13) = energy
164 ENDIF
165C
166 ivolu(14) = idef
167 rvolu(16) = vol
168 rvolu(12) = pres
169C
170 IF (ispmd+1==pmain) THEN
171 wfext=wfext+(half*(q+qold+pres+pold)-pext)*dv
172 fsav(1)=amtot
173 fsav(2)=vol
174 fsav(3)=pres
175 fsav(4)=area
176 fsav(5)=temperature
177 fsav(6)=zero
178 fsav(7)=zero
179 fsav(8)=zero
180 fsav(9)=zero
181 fsav(10)=zero
182 fsav(11)=zero
183 fsav(12)=gama
184 ENDIF
185C
186 rvolu(22) = zero
187 rvolu(23) = q
188 rvolu(24) = zero
189C
190 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine area(d1, x, x2, y, y2, eint, stif0)
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21

◆ volpvgb()

subroutine volpvgb ( integer, dimension(*) ivolu,
rvolu,
vol,
fsav,
integer nvent,
integer, dimension(nibhol,*) ibaghol,
rbaghol,
normal,
integer nn,
type (surf_), dimension(nsurf) igrsurf,
integer, dimension(nparg,*) iparg,
type(elbuf_struct_), dimension(ngroup), target elbuf_tab,
integer, dimension(*) fr_mv,
integer, dimension(*), intent(in) igroupc,
integer, dimension(*), intent(in) igrouptg )

Definition at line 205 of file volpvg.F.

208C-----------------------------------------------
209C MONITORED VOLUMES TYPE PERFECT GAS, INPUT STARTER >= 4.4 (FUITES).
210C-----------------------------------------------
211C M o d u l e s
212C-----------------------------------------------
213 USE elbufdef_mod
214 USE groupdef_mod
215C-----------------------------------------------
216C I m p l i c i t T y p e s
217C-----------------------------------------------
218#include "implicit_f.inc"
219C-----------------------------------------------
220C C o m m o n B l o c k s
221C-----------------------------------------------
222#include "param_c.inc"
223#include "units_c.inc"
224#include "com01_c.inc"
225#include "com04_c.inc"
226#include "com08_c.inc"
227#include "task_c.inc"
228#include "vect01_c.inc"
229C-----------------------------------------------
230C D u m m y A r g u m e n t s
231C-----------------------------------------------
232 INTEGER IVOLU(*),NVENT,IBAGHOL(NIBHOL,*), FR_MV(*),
233 . NN,IPARG(NPARG,*)
234 INTEGER, INTENT(IN) :: IGROUPC(*), IGROUPTG(*)
235 my_real
236 . rvolu(*),fsav(*),rbaghol(nrbhol,*),normal(3,*)
237 TYPE(ELBUF_STRUCT_), DIMENSION(NGROUP), TARGET :: ELBUF_TAB
238 TYPE (SURF_) , DIMENSION(NSURF) :: IGRSURF
239C-----------------------------------------------
240C L o c a l V a r i a b l e s
241C-----------------------------------------------
242 INTEGER K,IDEF,KK,PMAIN,
243 . II,IPVENT,NNC,KAD,IPORT,IPORP,IPORA,
244 . IOFF,NG,NEL,ISTRA,IEXPAN,IOK
245
246 my_real
247 . gama, pext, pdef, avent, vol, vinc,
248 . amtot, energy, p, ro,
249 . u, deout, dmout, tvent,area, pcrit, aoutot,
250 . apvent,aout,flout,de,deri,dtpdefi,dtpdefc,
251 . f1(nn),scalt,scalp,scals
252 my_real
253 . get_u_func
254 EXTERNAL get_u_func
255 DOUBLE PRECISION
256 . FRMV6(6)
257 my_real,
258 . DIMENSION(:), POINTER :: offg
259C-----------------------------------------------
260 pmain = fr_mv(nspmd+2)
261 idef = ivolu(14)
262 IF(idef>0)GOTO 999
263C--------------------------------
264 gama =rvolu(1)
265 pext =rvolu(3)
266 vinc =rvolu(5)
267 vol =rvolu(16)
268 energy =rvolu(13)
269 amtot =rvolu(20)
270 scalt =rvolu(26)
271 scalp =rvolu(27)
272 scals =rvolu(28)
273C
274 p =rvolu(12)
275 area =rvolu(18)
276C
277 ro = amtot/(vol-vinc)
278 pcrit = p*(two/(gama+one))**(gama/(gama-one))
279C
280C FLUX SORTANT PAR LES TROUS
281C
282 aoutot=zero
283 DO ii = 1, nvent
284 idef = ibaghol(1, ii)
285 ipvent = ibaghol(2, ii)
286C
287 pdef = rbaghol(1, ii)
288 dtpdefi= rbaghol(4, ii)
289 dtpdefc= rbaghol(5, ii)
290 avent = rbaghol(2, ii)
291 tvent = rbaghol(3, ii)
292C
293 IF(idef<=0.AND.p>pdef+pext.
294 . and.dtpdefc>dtpdefi.
295 . and.vol>0.001*area**1.5) THEN
296 idef=abs(idef)
297 IF(ispmd+1==pmain) THEN
298 WRITE(iout,*)
299 . ' *** MONITORED VOLUME MEMBRANE IS DEFLATED ***'
300 WRITE(iout,*)' *** MONITORED VOLUME ',ivolu(1),
301 . ' VENT HOLES MEMBRANE NUMBER ',ii,' ***'
302 WRITE(istdo,*)' *** VENT HOLES MEMBRANE IS DEFLATED ***'
303 ENDIF
304 ENDIF
305 IF(idef<=0 .AND. tt>tvent) THEN
306 idef=abs(idef)
307 IF(ispmd+1==pmain) THEN
308 WRITE(iout,*)
309 . ' *** MONITORED VOLUME VENTING STARTS ***'
310 WRITE(iout,*) ' *** MONITORED VOLUME ',ivolu(1),
311 . ' VENT HOLES MEMBRANE NUMBER ',ii,' ***'
312 WRITE(istdo,*)' *** VENTING STARTS ***'
313 ENDIF
314 ENDIF
315C
316 IF(ipvent/=0)THEN
317C
318 IF(idef==1)THEN
319c APVENT = ZERO
320 nnc=igrsurf(ipvent)%NSEG
321 DO kk=1,nnc
322 IF(igrsurf(ipvent)%ELTYP(kk)==3)THEN
323C segment from shell
324 k=igrsurf(ipvent)%ELEM(kk)
325 ELSEIF(igrsurf(ipvent)%ELTYP(kk)==7)THEN
326C segment from sh3n
327 k=igrsurf(ipvent)%ELEM(kk) + numelc
328 ELSE
329C segment only
330 k=igrsurf(ipvent)%ELEM(kk) + numelc + numeltg
331 ENDIF
332 f1(kk) = sqrt( normal(1,k)**2+normal(2,k)**2+normal(3,k)**2 )
333 kad=kad+nisx
334 ENDDO
335c AOUT=APVENT
336C
337 ELSEIF(idef>=2)THEN
338c APVENT = ZERO
339 nnc=igrsurf(ipvent)%NSEG
340 DO kk=1,nnc
341 IF(igrsurf(ipvent)%ELTYP(kk)==3)THEN
342C segment from shell
343 k=igrsurf(ipvent)%ELEM(kk)
344 iok = 0
345 ng=igroupc(k)
346 ity=iparg(5,ng)
347 IF(ity==3)THEN
348 iok = 1
349 ENDIF
350 IF (iok == 1) THEN
351 nel =iparg(2,ng)
352 nft =iparg(3,ng)
353 iad =iparg(4,ng)
354 npt =iparg(6,ng)
355 istra =iparg(44,ng)
356 jhbe =iparg(23,ng)
357 iexpan=iparg(49,ng)
358 offg => elbuf_tab(ng)%GBUF%OFF
359 ioff = int(offg(k-nft))
360 ELSE
361 ioff = 1
362 ENDIF
363 ELSEIF(igrsurf(ipvent)%ELTYP(kk)==7)THEN
364C segment from sh3n
365 k=igrsurf(ipvent)%ELEM(kk)
366 iok = 0
367 ng=igrouptg(k)
368 ity=iparg(5,ng)
369 IF(ity==7)THEN
370 iok = 1
371 ENDIF
372 IF (iok == 1) THEN
373 nel =iparg(2,ng)
374 nft =iparg(3,ng)
375 iad =iparg(4,ng)
376 npt =iparg(6,ng)
377 istra =iparg(44,ng)
378 jhbe =iparg(23,ng)
379 iexpan=iparg(49,ng)
380 offg => elbuf_tab(ng)%GBUF%OFF
381 ioff=int(offg(k-nft))
382 ELSE
383 ioff = 1
384 ENDIF
385 k=k+numelc
386 ELSE
387C segment only
388 ioff=1
389 ENDIF
390 IF(ioff==0) THEN
391 f1(kk) = sqrt( normal(1,k)**2+normal(2,k)**2+normal(3,k)**2 )
392 ELSE
393 f1(kk) = zero
394 END IF
395 kad=kad+nisx
396 ENDDO
397c AOUT=APVENT
398 ENDIF
399C
400
401 IF (idef==1.OR.idef>=2) THEN
402 DO kk = 1, 6
403 frmv6(kk) = zero
404 END DO
405 CALL sum_6_float(1, nnc, f1, frmv6, 1)
406 IF (nspmd>1) THEN
407 IF(fr_mv(ispmd+1)/=0) THEN
408 CALL spmd_exch_fr6(fr_mv,frmv6,6)
409 END IF
410 ENDIF
411 apvent = frmv6(1)+frmv6(2)+frmv6(3)+
412 . frmv6(4)+frmv6(5)+frmv6(6)
413C
414 aout = apvent
415 ELSE
416 aout =avent
417 avent=1.0
418 END IF
419 ELSE
420 aout =avent
421 avent=1.0
422 ENDIF
423C
424 IF(idef>0 .AND. p>pext.
425 . and.vol>em3*area**1.5)THEN
426 iport =ibaghol(3,ii)
427 iporp =ibaghol(4,ii)
428 ipora =ibaghol(5,ii)
429 IF(ipora/=0.AND.ipvent/=0)THEN
430 aout=avent*get_u_func(ipora,aout*scals,deri)
431 ELSE
432 aout=avent*aout
433 ENDIF
434 IF(iport/=0)aout=aout*get_u_func(iport,tt*scalt,deri)
435 IF(iporp/=0)aout=aout*get_u_func(iporp,(p-pext)*scalp,deri)
436 aoutot=aoutot+aout
437 ENDIF
438 ibaghol(1,ii)=idef
439 ENDDO
440C-------
441 IF(aoutot>0.)THEN
442 pext = max(pext,pcrit)
443 u=two*gama/(gama-one)*p/ro*(one-(pext/p)**((gama-one)/gama))
444 u=sqrt(u)
445 de=(energy/(vol-vinc)+p)*(pext/p)**(one/gama)
446 u=min(u,(p-pext)*half*(vol-vinc)
447 . /(gama-one)/de/max(em20,aoutot*dt1))
448 u=min(u,half*(vol-vinc)/max(em20,aoutot*dt1))
449 flout=aoutot*u
450 deout=flout*de
451 dmout=flout*ro*(pext/p)**(one/gama)
452 ELSE
453 u=zero
454 deout=zero
455 dmout=zero
456 flout= zero
457 ENDIF
458C-------
459 IF(ispmd+1==pmain) THEN
460 fsav(6)=aoutot
461 fsav(7)=flout/max(em20,aoutot)
462 ENDIF
463C
464 rvolu(22)=rvolu(22) + deout
465 rvolu(24)=rvolu(24) + dmout
466C--------------------------------
467 999 CONTINUE
468 RETURN
subroutine sum_6_float(jft, jlt, f, f6, n)
Definition parit.F:64
subroutine spmd_exch_fr6(fr, fs6, len)