xref: /phasta/phSolver/compressible/e3ql.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine e3ql (ycl,      shp,     shgl,
2*59599516SKenneth E. Jansen     &                   xl,      ql,      xmudmi,
3*59599516SKenneth E. Jansen     &                   sgn )
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansenc----------------------------------------------------------------------
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc This routine computes the local diffusive flux vector using a
8*59599516SKenneth E. Jansenc local projection algorithm
9*59599516SKenneth E. Jansenc
10*59599516SKenneth E. Jansenc input:
11*59599516SKenneth E. Jansenc  ycl     (npro,nshl,ndof)      : Y variables
12*59599516SKenneth E. Jansenc  shp    (nen,ngauss)           : element shape-functions
13*59599516SKenneth E. Jansenc  shgl   (nsd,nen,ngauss)       : element local-grad-shape-functions
14*59599516SKenneth E. Jansenc  xl     (npro,nshape,nsd)      : nodal coordinates at current step
15*59599516SKenneth E. Jansenc  sgn    (npro,nshl)            : signs for reversed shape functions
16*59599516SKenneth E. Jansenc
17*59599516SKenneth E. Jansenc output:
18*59599516SKenneth E. Jansenc  ql     (npro,nshape,(nflow-1)*nsd) : element RHS diffusion residual
19*59599516SKenneth E. Jansenc
20*59599516SKenneth E. Jansenc----------------------------------------------------------------------
21*59599516SKenneth E. Jansenc
22*59599516SKenneth E. Jansen        include "common.h"
23*59599516SKenneth E. Jansenc
24*59599516SKenneth E. Jansen        dimension ycl(npro,nshl,ndof),
25*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),
26*59599516SKenneth E. Jansen     &            shgl(nsd,nshl,ngauss),
27*59599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),       sgn(npro,nshl),
28*59599516SKenneth E. Jansen     &            ql(npro,nshl,idflx), xmudmi(npro,ngauss)
29*59599516SKenneth E. Jansenc
30*59599516SKenneth E. Jansenc local arrays
31*59599516SKenneth E. Jansenc
32*59599516SKenneth E. Jansen        dimension g1yi(npro,nflow),          g2yi(npro,nflow),
33*59599516SKenneth E. Jansen     &            g3yi(npro,nflow),          shg(npro,nshl,nsd),
34*59599516SKenneth E. Jansen     &            dxidx(npro,nsd,nsd),       WdetJ(npro),
35*59599516SKenneth E. Jansen     &            T(npro),                   cp(npro),
36*59599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
37*59599516SKenneth E. Jansen     &            u3(npro),                  rmu(npro),
38*59599516SKenneth E. Jansen     &            rlm(npro),                 rlm2mu(npro),
39*59599516SKenneth E. Jansen     &            con(npro),                 rminv(npro,nshl,nshl),
40*59599516SKenneth E. Jansen     &            qrl(npro,nshl,(nflow-1)*nsd)
41*59599516SKenneth E. Jansenc
42*59599516SKenneth E. Jansen        dimension qdi(npro,nsd*(nflow-1)),   shape(npro,nshl),
43*59599516SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),      indx(nshl),
44*59599516SKenneth E. Jansen     &            rmass(npro,nshl,nshl),     rho(npro)
45*59599516SKenneth E. Jansen
46*59599516SKenneth E. Jansen        real*8 tmp(npro)
47*59599516SKenneth E. Jansenc
48*59599516SKenneth E. Jansenc.... loop through the integration points
49*59599516SKenneth E. Jansenc
50*59599516SKenneth E. Jansen        rminv = zero
51*59599516SKenneth E. Jansen        rmass = zero
52*59599516SKenneth E. Jansen        qrl = zero
53*59599516SKenneth E. Jansen
54*59599516SKenneth E. Jansen        do intp = 1, ngauss
55*59599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle ! precaution
56*59599516SKenneth E. Jansenc
57*59599516SKenneth E. Jansenc.... create a matrix of shape functions (and derivatives) for each
58*59599516SKenneth E. Jansenc     element at this quadrature point. These arrays will contain
59*59599516SKenneth E. Jansenc     the correct signs for the hierarchic basis
60*59599516SKenneth E. Jansenc
61*59599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
62*59599516SKenneth E. Jansen     &              shape,        shdrv)
63*59599516SKenneth E. Jansenc
64*59599516SKenneth E. Jansenc
65*59599516SKenneth E. Jansenc.... initialize
66*59599516SKenneth E. Jansenc
67*59599516SKenneth E. Jansen        qdi = zero
68*59599516SKenneth E. Jansenc
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansenc.... calculate the integration variables necessary for the
71*59599516SKenneth E. Jansenc     formation of q
72*59599516SKenneth E. Jansenc
73*59599516SKenneth E. Jansen
74*59599516SKenneth E. Jansen        call e3qvar   (ycl,       shape,        shdrv,
75*59599516SKenneth E. Jansen     &                 rho,       xl,           g1yi,
76*59599516SKenneth E. Jansen     &                 g2yi,      g3yi,         shg,
77*59599516SKenneth E. Jansen     &                 dxidx,     WdetJ,        T,
78*59599516SKenneth E. Jansen     &                 cp,        u1,           u2,
79*59599516SKenneth E. Jansen     &                 u3                                 )
80*59599516SKenneth E. Jansenc
81*59599516SKenneth E. Jansenc
82*59599516SKenneth E. Jansenc.... compute diffusive flux vector at this integration point
83*59599516SKenneth E. Jansenc
84*59599516SKenneth E. Jansenc
85*59599516SKenneth E. Jansenc.... get material properties
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansen        call getDiff (T,        cp,    rho,      ycl,
88*59599516SKenneth E. Jansen     &               rmu,      rlm,    rlm2mu,   con, shape,
89*59599516SKenneth E. Jansen     &               xmudmi,   xl)
90*59599516SKenneth E. Jansenc
91*59599516SKenneth E. Jansenc.... compute diffusive fluxes
92*59599516SKenneth E. Jansenc
93*59599516SKenneth E. Jansenc.... diffusive flux in x1-direction
94*59599516SKenneth E. Jansenc
95*59599516SKenneth E. Jansen        qdi(:,1) =  rlm2mu      * g1yi(:,2)
96*59599516SKenneth E. Jansen     &               +      rlm * g2yi(:,3)
97*59599516SKenneth E. Jansen     &               +      rlm * g3yi(:,4)
98*59599516SKenneth E. Jansen        qdi(:,2) =  rmu         * g1yi(:,3)
99*59599516SKenneth E. Jansen     &               +      rmu * g2yi(:,2)
100*59599516SKenneth E. Jansen        qdi(:,3) =  rmu         * g1yi(:,4)
101*59599516SKenneth E. Jansen     &               +      rmu * g3yi(:,2)
102*59599516SKenneth E. Jansen        qdi(:,4) =  rlm2mu * u1 * g1yi(:,2) + rmu * u2 * g1yi(:,3)
103*59599516SKenneth E. Jansen     &                                   +    rmu * u3 * g1yi(:,4)
104*59599516SKenneth E. Jansen     &               +    rmu * u2 * g2yi(:,2) + rlm * u1 * g2yi(:,3)
105*59599516SKenneth E. Jansen     &               +    rmu * u3 * g3yi(:,2) + rlm * u1 * g3yi(:,4)
106*59599516SKenneth E. Jansen     &               +    con      * g1yi(:,5)
107*59599516SKenneth E. Jansen
108*59599516SKenneth E. Jansenc
109*59599516SKenneth E. Jansenc.... diffusive flux in x2-direction
110*59599516SKenneth E. Jansenc
111*59599516SKenneth E. Jansen        qdi(:,5) =       rmu * g1yi(:,3)
112*59599516SKenneth E. Jansen     &            +      rmu * g2yi(:,2)
113*59599516SKenneth E. Jansen        qdi(:,6) =       rlm * g1yi(:,2)
114*59599516SKenneth E. Jansen     &            +   rlm2mu * g2yi(:,3)
115*59599516SKenneth E. Jansen     &            +      rlm * g3yi(:,4)
116*59599516SKenneth E. Jansen        qdi(:,7) =       rmu * g2yi(:,4)
117*59599516SKenneth E. Jansen     &            +      rmu * g3yi(:,3)
118*59599516SKenneth E. Jansen        qdi(:,8) =  rlm * u2 * g1yi(:,2) +    rmu * u1 * g1yi(:,3)
119*59599516SKenneth E. Jansen     &            + rmu * u1 * g2yi(:,2) + rlm2mu * u2 * g2yi(:,3)
120*59599516SKenneth E. Jansen     &            + rmu * u3 * g2yi(:,4)
121*59599516SKenneth E. Jansen     &            + rmu * u3 * g3yi(:,3) +    rlm * u2 * g3yi(:,4)
122*59599516SKenneth E. Jansen     &            +    con      * g2yi(:,5)
123*59599516SKenneth E. Jansenc
124*59599516SKenneth E. Jansenc.... diffusive flux in x3-direction
125*59599516SKenneth E. Jansenc
126*59599516SKenneth E. Jansen        qdi(:,9 ) =       rmu * g1yi(:,4)
127*59599516SKenneth E. Jansen     &            +      rmu * g3yi(:,2)
128*59599516SKenneth E. Jansen        qdi(:,10) =       rmu * g2yi(:,4)
129*59599516SKenneth E. Jansen     &            +      rmu * g3yi(:,3)
130*59599516SKenneth E. Jansen        qdi(:,11) =       rlm * g1yi(:,2)
131*59599516SKenneth E. Jansen     &            +      rlm * g2yi(:,3)
132*59599516SKenneth E. Jansen     &            +   rlm2mu * g3yi(:,4)
133*59599516SKenneth E. Jansen        qdi(:,12) =     rlm * u3 * g1yi(:,2) + rmu * u1 * g1yi(:,4)
134*59599516SKenneth E. Jansen     &            +    rlm * u3 * g2yi(:,3) + rmu * u2 * g2yi(:,4)
135*59599516SKenneth E. Jansen     &            +    rmu * u1 * g3yi(:,2) + rmu * u2 * g3yi(:,3)
136*59599516SKenneth E. Jansen     &            + rlm2mu * u3 * g3yi(:,4)
137*59599516SKenneth E. Jansen     &            +    con      * g3yi(:,5)
138*59599516SKenneth E. Jansen
139*59599516SKenneth E. Jansenc
140*59599516SKenneth E. Jansenc
141*59599516SKenneth E. Jansenc.... assemble contribution of qdi to qrl,i.e., contribution to
142*59599516SKenneth E. Jansenc     each element shape function
143*59599516SKenneth E. Jansenc
144*59599516SKenneth E. Jansen        tmp = Qwt(lcsyst,intp)
145*59599516SKenneth E. Jansen        if ( lcsyst .eq. 1) then
146*59599516SKenneth E. Jansen     	   tmp = tmp*(three/four)
147*59599516SKenneth E. Jansen        end if
148*59599516SKenneth E. Jansen
149*59599516SKenneth E. Jansen        do i=1,nshl
150*59599516SKenneth E. Jansen           qrl(:,i,1 ) = qrl(:,i,1 )+ shape(:,i)*tmp*qdi(:,1 )
151*59599516SKenneth E. Jansen           qrl(:,i,2 ) = qrl(:,i,2 )+ shape(:,i)*tmp*qdi(:,2 )
152*59599516SKenneth E. Jansen           qrl(:,i,3 ) = qrl(:,i,3 )+ shape(:,i)*tmp*qdi(:,3 )
153*59599516SKenneth E. Jansen           qrl(:,i,4 ) = qrl(:,i,4 )+ shape(:,i)*tmp*qdi(:,4 )
154*59599516SKenneth E. Jansen           qrl(:,i,5 ) = qrl(:,i,5 )+ shape(:,i)*tmp*qdi(:,5 )
155*59599516SKenneth E. Jansen           qrl(:,i,6 ) = qrl(:,i,6 )+ shape(:,i)*tmp*qdi(:,6 )
156*59599516SKenneth E. Jansen           qrl(:,i,7 ) = qrl(:,i,7 )+ shape(:,i)*tmp*qdi(:,7 )
157*59599516SKenneth E. Jansen           qrl(:,i,8 ) = qrl(:,i,8 )+ shape(:,i)*tmp*qdi(:,8 )
158*59599516SKenneth E. Jansen           qrl(:,i,9 ) = qrl(:,i,9 )+ shape(:,i)*tmp*qdi(:,9 )
159*59599516SKenneth E. Jansen           qrl(:,i,10) = qrl(:,i,10)+ shape(:,i)*tmp*qdi(:,10)
160*59599516SKenneth E. Jansen           qrl(:,i,11) = qrl(:,i,11)+ shape(:,i)*tmp*qdi(:,11)
161*59599516SKenneth E. Jansen           qrl(:,i,12) = qrl(:,i,12)+ shape(:,i)*tmp*qdi(:,12)
162*59599516SKenneth E. Jansen        enddo
163*59599516SKenneth E. Jansenc
164*59599516SKenneth E. Jansenc.... add contribution to local mass matrix
165*59599516SKenneth E. Jansenc
166*59599516SKenneth E. Jansen        do i=1,nshl
167*59599516SKenneth E. Jansen           do j=1,nshl
168*59599516SKenneth E. Jansen              rmass(:,i,j) = rmass(:,i,j)+shape(:,i)*shape(:,j)*tmp
169*59599516SKenneth E. Jansen           enddo
170*59599516SKenneth E. Jansen        enddo
171*59599516SKenneth E. Jansenc
172*59599516SKenneth E. Jansenc.... end of the loop over integration points
173*59599516SKenneth E. Jansenc
174*59599516SKenneth E. Jansen       enddo
175*59599516SKenneth E. Jansen       if ( lcsyst .eq. 1) then
176*59599516SKenneth E. Jansen          qrl = qrl/6.d0
177*59599516SKenneth E. Jansen          rmass = rmass/6.0
178*59599516SKenneth E. Jansen       end if
179*59599516SKenneth E. Jansenc
180*59599516SKenneth E. Jansenc.... find the inverse of the local mass matrix for each element
181*59599516SKenneth E. Jansenc
182*59599516SKenneth E. Jansen       do iel=1,npro
183*59599516SKenneth E. Jansen
184*59599516SKenneth E. Jansen          do i=1,nshl   ! form the identy matrix
185*59599516SKenneth E. Jansen             do j=1,nshl
186*59599516SKenneth E. Jansen                rminv(iel,i,j) = 0.0
187*59599516SKenneth E. Jansen             enddo
188*59599516SKenneth E. Jansen             rminv(iel,i,i)=1.0
189*59599516SKenneth E. Jansen          enddo
190*59599516SKenneth E. Jansenc
191*59599516SKenneth E. Jansenc.... LU factor the mass matrix
192*59599516SKenneth E. Jansenc
193*59599516SKenneth E. Jansen          call ludcmp(rmass(iel,:,:),nshl,nshl,indx,d)
194*59599516SKenneth E. Jansenc
195*59599516SKenneth E. Jansenc.... back substitute with the identy matrix to find the
196*59599516SKenneth E. Jansenc     matrix inverse
197*59599516SKenneth E. Jansenc
198*59599516SKenneth E. Jansen          do j=1,nshl
199*59599516SKenneth E. Jansen             call lubksb(rmass(iel,:,:),nshl,nshl,
200*59599516SKenneth E. Jansen     &                                  indx,rminv(iel,:,j))
201*59599516SKenneth E. Jansen          enddo
202*59599516SKenneth E. Jansen       enddo
203*59599516SKenneth E. Jansen
204*59599516SKenneth E. Jansenc
205*59599516SKenneth E. Jansenc.... find the modal coefficients of ql by multiplying by the inverse of
206*59599516SKenneth E. Jansenc     the local mass matrix
207*59599516SKenneth E. Jansenc
208*59599516SKenneth E. Jansen       do iel=1,npro
209*59599516SKenneth E. Jansen          do j=1,12
210*59599516SKenneth E. Jansen             ql(iel,:,j) = matmul( rminv(iel,:,:),qrl(iel,:,j) )
211*59599516SKenneth E. Jansen          enddo
212*59599516SKenneth E. Jansen       enddo
213*59599516SKenneth E. Jansenc
214*59599516SKenneth E. Jansenc.... return
215*59599516SKenneth E. Jansenc
216*59599516SKenneth E. Jansen       return
217*59599516SKenneth E. Jansen       end
218