CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
calcrhofacomp.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine calcrhofacomp (nface, vfa, shcon, ielmat, ntmat_, mi, ielfa, ipnei, vel, nef, flux, gradpel, gradtel, xxj, betam, xlet)
 

Function/Subroutine Documentation

◆ calcrhofacomp()

subroutine calcrhofacomp ( integer  nface,
real*8, dimension(0:7,*)  vfa,
real*8, dimension(0:3,ntmat_,*)  shcon,
integer, dimension(mi(3),*)  ielmat,
integer  ntmat_,
integer, dimension(*)  mi,
integer, dimension(4,*)  ielfa,
integer, dimension(*)  ipnei,
real*8, dimension(nef,0:7)  vel,
integer  nef,
real*8, dimension(*)  flux,
real*8, dimension(3,*)  gradpel,
real*8, dimension(3,*)  gradtel,
real*8, dimension(3,*)  xxj,
real*8  betam,
real*8, dimension(*)  xlet 
)
22 !
23 ! calculation of the density at the face centers
24 ! (compressible fluids)
25 !
26  implicit none
27 !
28  integer nface,i,j,imat,ntmat_,mi(*),ipnei(*),nef,iel1,iel2,
29  & ielmat(mi(3),*),ielfa(4,*),indexf
30 !
31  real*8 t1l,vfa(0:7,*),shcon(0:3,ntmat_,*),vel(nef,0:7),flux(*),
32  & r,gradpel(3,*),gradtel(3,*),xxj(3,*),gamma,betam,phic,vud,
33  & vcd,xlet(*)
34 !
35 c$omp parallel default(none)
36 c$omp& shared(nface,vfa,ielmat,ielfa,shcon,ipnei,flux,gradpel,gradtel,
37 c$omp& xxj,betam,xlet,vel)
38 c$omp& private(i,j,t1l,imat,iel1,iel2,indexf,vcd,vud,r,gamma,phic)
39 c$omp do
40  do i=1,nface
41  t1l=vfa(0,i)
42 !
43 ! take the material of the first adjacent element
44 !
45  imat=ielmat(1,ielfa(1,i))
46  r=shcon(3,1,imat)
47 !
48 ! specific gas constant
49 !
50  vfa(5,i)=vfa(4,i)/(r*t1l)
51 !
52 ! calculate gamma
53 !
54  iel2=ielfa(2,i)
55 !
56 ! faces with only one neighbor need not be treated
57 !
58  if(iel2.le.0) cycle
59  iel1=ielfa(1,i)
60  j=ielfa(4,i)
61  indexf=ipnei(iel1)+j
62 !
63  vcd=vel(iel2,5)-vel(iel1,5)
64  if(dabs(vcd).lt.1.d-3*dabs(vel(iel1,5))) vcd=0.d0
65 !
66  if(flux(indexf).ge.0.d0) then
67  vud=2.d0*xlet(indexf)*
68  & ((gradpel(1,iel1)-vel(iel1,5)*r*gradtel(1,iel1))
69  & *xxj(1,indexf)+
70  & (gradpel(2,iel1)-vel(iel1,5)*r*gradtel(2,iel1))
71  & *xxj(2,indexf)+
72  & (gradpel(3,iel1)-vel(iel1,5)*r*gradtel(3,iel1))
73  & *xxj(3,indexf))/(r*vel(iel1,0))
74  else
75  vud=2.d0*xlet(indexf)*
76  & ((gradpel(1,iel2)-vel(iel2,5)*r*gradtel(1,iel2))
77  & *xxj(1,indexf)+
78  & (gradpel(2,iel2)-vel(iel2,5)*r*gradtel(2,iel2))
79  & *xxj(2,indexf)+
80  & (gradpel(3,iel2)-vel(iel2,5)*r*gradtel(3,iel2))
81  & *xxj(3,indexf))/(r*vel(iel2,0))
82  endif
83 !
84  if(dabs(vud).lt.1.d-20) then
85  gamma=0.d0
86  cycle
87  endif
88 !
89  phic=1.d0-vcd/vud
90 !
91  if(phic.ge.1.d0) then
92  gamma=0.d0
93  elseif(phic.le.0.d0) then
94  gamma=0.d0
95  elseif(betam.le.phic) then
96  gamma=1.d0
97  else
98  gamma=phic/betam
99  endif
100 !
101 ! mixture of central difference and upwind difference
102 !
103  if(flux(indexf).ge.0.d0) then
104  vfa(5,i)=gamma*vfa(5,i)+(1.d0-gamma)*vel(iel1,5)
105  else
106  vfa(5,i)=gamma*vfa(5,i)+(1.d0-gamma)*vel(iel2,5)
107  endif
108 !
109 c write(*,*) 'calcrhofacomp ',i,gamma
110  enddo
111 c$omp end do
112 c$omp end parallel
113 !
114  return
subroutine flux(node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, kflag, v, xflow, f, nodef, idirf, df, cp, R, rho, physcon, g, co, dvi, numf, vold, set, shcon, nshcon, rhcon, nrhcon, ntmat_, mi, ider, ttime, time, iaxial)
Definition: flux.f:24
Hosted by OpenAircraft.com, (Michigan UAV, LLC)