xref: /phasta/phSolver/compressible/e3qvar.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32) !
1*59599516SKenneth E. Jansen        subroutine e3qvar (ycl, shp,     shgl,    rho,
2*59599516SKenneth E. Jansen     &                     xl,  g1yi,    g2yi,    g3yi,
3*59599516SKenneth E. Jansen     &                     shg, dxidx,   WdetJ,   T,
4*59599516SKenneth E. Jansen     &                     cp,  u1,      u2,      u3)
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc----------------------------------------------------------------------
7*59599516SKenneth E. Jansenc
8*59599516SKenneth E. Jansenc  This routine computes the variables at integration point
9*59599516SKenneth E. Jansenc  necessary for the computation of the diffusive flux vector.
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansenc input:
12*59599516SKenneth E. Jansenc  ycl    (npro,nshape,ndof)      : primitive variables
13*59599516SKenneth E. Jansenc  shp    (npro,nshape)         : element shape-functions
14*59599516SKenneth E. Jansenc  shgl   (npro,nsd,nshape)     : element local-grad-shape-functions
15*59599516SKenneth E. Jansenc  xl     (npro,nenl,nsd)       : nodal coordinates at current step
16*59599516SKenneth E. Jansenc
17*59599516SKenneth E. Jansenc output:
18*59599516SKenneth E. Jansenc  g1yi   (npro,nflow)           : grad-y in direction 1
19*59599516SKenneth E. Jansenc  g2yi   (npro,nflow)           : grad-y in direction 2
20*59599516SKenneth E. Jansenc  g3yi   (npro,nflow)           : grad-y in direction 3
21*59599516SKenneth E. Jansenc  shg    (npro,nshape,nsd)       : element global grad-shape-functions
22*59599516SKenneth E. Jansenc  dxidx  (npro,nsd,nsd)        : inverse of deformation gradient
23*59599516SKenneth E. Jansenc  WdetJ  (npro)                : weighted Jacobian
24*59599516SKenneth E. Jansenc  T      (npro)                : temperature
25*59599516SKenneth E. Jansenc  cp     (npro)                : specific heat at constant pressure
26*59599516SKenneth E. Jansenc  u1     (npro)                : x1-velocity component
27*59599516SKenneth E. Jansenc  u2     (npro)                : x2-velocity component
28*59599516SKenneth E. Jansenc  u3     (npro)                : x3-velocity component
29*59599516SKenneth E. Jansenc
30*59599516SKenneth E. Jansenc
31*59599516SKenneth E. Jansenc Zdenek Johan, Summer 1990. (Modified from e2ivar.f)
32*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991. (Fortran 90)
33*59599516SKenneth E. Jansenc Kenneth Jansen, Winter 1997. Primitive Variables
34*59599516SKenneth E. Jansenc----------------------------------------------------------------------
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansen        include "common.h"
37*59599516SKenneth E. Jansenc
38*59599516SKenneth E. Jansenc  passed arrays
39*59599516SKenneth E. Jansenc
40*59599516SKenneth E. Jansen        dimension ycl(npro,nshl,ndof),  rho(npro),
41*59599516SKenneth E. Jansen     &            shp(npro,nshl),      tmp(npro),
42*59599516SKenneth E. Jansen     &            shgl(npro,nsd,nshl), xl(npro,nenl,nsd),
43*59599516SKenneth E. Jansen     &            g1yi(npro,nflow),    g2yi(npro,nflow),
44*59599516SKenneth E. Jansen     &            g3yi(npro,nflow),    shg(npro,nshl,nsd),
45*59599516SKenneth E. Jansen     &            dxidx(npro,nsd,nsd), WdetJ(npro),
46*59599516SKenneth E. Jansen     &            T(npro),             cp(npro),
47*59599516SKenneth E. Jansen     &            u1(npro),            u2(npro),
48*59599516SKenneth E. Jansen     &            u3(npro),            pres(npro)
49*59599516SKenneth E. Jansenc
50*59599516SKenneth E. Jansenc  local arrays
51*59599516SKenneth E. Jansenc
52*59599516SKenneth E. Jansen        dimension dxdxi(npro,nsd,nsd), Sclr(npro)
53*59599516SKenneth E. Jansen
54*59599516SKenneth E. Jansen        ttim(34) = ttim(34) - secs(0.0)
55*59599516SKenneth E. Jansenc
56*59599516SKenneth E. Jansenc.... ------------->  Primitive variables at int. point  <--------------
57*59599516SKenneth E. Jansenc
58*59599516SKenneth E. Jansenc.... compute primitive variables
59*59599516SKenneth E. Jansenc
60*59599516SKenneth E. Jansen       pres = zero
61*59599516SKenneth E. Jansen       u1   = zero
62*59599516SKenneth E. Jansen       u2   = zero
63*59599516SKenneth E. Jansen       u3   = zero
64*59599516SKenneth E. Jansen       T    = zero
65*59599516SKenneth E. Jansen
66*59599516SKenneth E. Jansen
67*59599516SKenneth E. Jansenc
68*59599516SKenneth E. Jansen       do n = 1, nshl
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansenc  y(intp)=SUM_{a=1}^nshape (N_a(intp) Ya) (we don't need pressure here)
71*59599516SKenneth E. Jansenc
72*59599516SKenneth E. Jansen          pres = pres + shp(:,n) * ycl(:,n,1)
73*59599516SKenneth E. Jansen          u1   = u1   + shp(:,n) * ycl(:,n,2)
74*59599516SKenneth E. Jansen          u2   = u2   + shp(:,n) * ycl(:,n,3)
75*59599516SKenneth E. Jansen          u3   = u3   + shp(:,n) * ycl(:,n,4)
76*59599516SKenneth E. Jansen          T    = T    + shp(:,n) * ycl(:,n,5)
77*59599516SKenneth E. Jansen       enddo
78*59599516SKenneth E. Jansenc
79*59599516SKenneth E. Jansenc$$$c.... compute cp, only thermodynamic property needed for qdiff
80*59599516SKenneth E. Jansenc$$$c
81*59599516SKenneth E. Jansenc$$$       cp    = Rgas * gamma / gamma1
82*59599516SKenneth E. Jansenc.... get the thermodynamic properties
83*59599516SKenneth E. Jansenc
84*59599516SKenneth E. Jansen
85*59599516SKenneth E. Jansen        if (iLSet .ne. 0)then
86*59599516SKenneth E. Jansen           Sclr = zero
87*59599516SKenneth E. Jansen           isc=abs(iRANS)+6
88*59599516SKenneth E. Jansen           do n = 1, nshl
89*59599516SKenneth E. Jansen              Sclr = Sclr + shp(:,n) * ycl(:,n,isc)
90*59599516SKenneth E. Jansen           enddo
91*59599516SKenneth E. Jansen        endif
92*59599516SKenneth E. Jansenc
93*59599516SKenneth E. Jansen        if (Navier .eq. 1) ithm = 7
94*59599516SKenneth E. Jansen
95*59599516SKenneth E. Jansen        call getthm (pres,            T,                  Sclr,
96*59599516SKenneth E. Jansen     &               tmp,             rho,                tmp,
97*59599516SKenneth E. Jansen     &               tmp,             tmp,                tmp,
98*59599516SKenneth E. Jansen     &               cp,              tmp,                tmp,
99*59599516SKenneth E. Jansen     &               tmp,             tmp)
100*59599516SKenneth E. Jansenc
101*59599516SKenneth E. Jansenc.... --------------------->  Element Metrics  <-----------------------
102*59599516SKenneth E. Jansenc
103*59599516SKenneth E. Jansen        call e3metric( xl,         shgl,        dxidx,
104*59599516SKenneth E. Jansen     &                 shg,        WdetJ)
105*59599516SKenneth E. Jansen
106*59599516SKenneth E. Jansenc
107*59599516SKenneth E. Jansenc.... --------------------->  Global Gradients  <-----------------------
108*59599516SKenneth E. Jansenc
109*59599516SKenneth E. Jansen        g1yi = zero
110*59599516SKenneth E. Jansen        g2yi = zero
111*59599516SKenneth E. Jansen        g3yi = zero
112*59599516SKenneth E. Jansenc
113*59599516SKenneth E. Jansenc
114*59599516SKenneth E. Jansen        do n = 1, nshl
115*59599516SKenneth E. Jansenc.... compute the global gradient of Y-variables
116*59599516SKenneth E. Jansenc
117*59599516SKenneth E. Jansenc
118*59599516SKenneth E. Jansenc  Y_{,x_i}=SUM_{a=1}^nenl (N_{a,x_i}(intp) Ya)
119*59599516SKenneth E. Jansenc
120*59599516SKenneth E. Jansen          g1yi(:,1) = g1yi(:,1) + shg(:,n,1) * ycl(:,n,1)
121*59599516SKenneth E. Jansen          g1yi(:,2) = g1yi(:,2) + shg(:,n,1) * ycl(:,n,2)
122*59599516SKenneth E. Jansen          g1yi(:,3) = g1yi(:,3) + shg(:,n,1) * ycl(:,n,3)
123*59599516SKenneth E. Jansen          g1yi(:,4) = g1yi(:,4) + shg(:,n,1) * ycl(:,n,4)
124*59599516SKenneth E. Jansen          g1yi(:,5) = g1yi(:,5) + shg(:,n,1) * ycl(:,n,5)
125*59599516SKenneth E. Jansenc
126*59599516SKenneth E. Jansen          g2yi(:,1) = g2yi(:,1) + shg(:,n,2) * ycl(:,n,1)
127*59599516SKenneth E. Jansen          g2yi(:,2) = g2yi(:,2) + shg(:,n,2) * ycl(:,n,2)
128*59599516SKenneth E. Jansen          g2yi(:,3) = g2yi(:,3) + shg(:,n,2) * ycl(:,n,3)
129*59599516SKenneth E. Jansen          g2yi(:,4) = g2yi(:,4) + shg(:,n,2) * ycl(:,n,4)
130*59599516SKenneth E. Jansen          g2yi(:,5) = g2yi(:,5) + shg(:,n,2) * ycl(:,n,5)
131*59599516SKenneth E. Jansenc
132*59599516SKenneth E. Jansen          g3yi(:,1) = g3yi(:,1) + shg(:,n,3) * ycl(:,n,1)
133*59599516SKenneth E. Jansen          g3yi(:,2) = g3yi(:,2) + shg(:,n,3) * ycl(:,n,2)
134*59599516SKenneth E. Jansen          g3yi(:,3) = g3yi(:,3) + shg(:,n,3) * ycl(:,n,3)
135*59599516SKenneth E. Jansen          g3yi(:,4) = g3yi(:,4) + shg(:,n,3) * ycl(:,n,4)
136*59599516SKenneth E. Jansen          g3yi(:,5) = g3yi(:,5) + shg(:,n,3) * ycl(:,n,5)
137*59599516SKenneth E. Jansen
138*59599516SKenneth E. Jansen        enddo
139*59599516SKenneth E. Jansen
140*59599516SKenneth E. Jansen        ttim(34) = ttim(34) + secs(0.0)
141*59599516SKenneth E. Jansen
142*59599516SKenneth E. Jansenc
143*59599516SKenneth E. Jansenc.... return
144*59599516SKenneth E. Jansenc
145*59599516SKenneth E. Jansen
146*59599516SKenneth E. Jansen       return
147*59599516SKenneth E. Jansen       end
148