xref: /phasta/phSolver/incompressible/e3ql.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine e3ql (yl,      dwl,     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  yl     (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,nshl,nsd*nsd) : element RHS diffusion residual
19*59599516SKenneth E. Jansenc
20*59599516SKenneth E. Jansenc----------------------------------------------------------------------
21*59599516SKenneth E. Jansenc
22*59599516SKenneth E. Jansen      use local_mass
23*59599516SKenneth E. Jansen      include "common.h"
24*59599516SKenneth E. Jansenc
25*59599516SKenneth E. Jansen      dimension yl(npro,nshl,ndof),        dwl(npro,nshl),
26*59599516SKenneth E. Jansen     &          shp(nshl,ngauss),          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,ndof),           g2yi(npro,ndof),
33*59599516SKenneth E. Jansen     &          g3yi(npro,ndof),           shg(npro,nshl,nsd),
34*59599516SKenneth E. Jansen     &          dxidx(npro,nsd,nsd),       WdetJ(npro),
35*59599516SKenneth E. Jansen     &          rmu(npro),
36*59599516SKenneth E. Jansen     &          rminv(npro,nshl,nshl),
37*59599516SKenneth E. Jansen     &          qrl(npro,nshl,nsd*nsd)
38*59599516SKenneth E. Jansenc
39*59599516SKenneth E. Jansen      dimension qdi(npro,nsd*nsd),    shape(npro,nshl),
40*59599516SKenneth E. Jansen     &          shdrv(npro,nsd,nshl),      indx(nshl),
41*59599516SKenneth E. Jansen     &          rmass(npro,nshl,nshl)
42*59599516SKenneth E. Jansen
43*59599516SKenneth E. Jansen
44*59599516SKenneth E. Jansen        real*8 tmp(npro)
45*59599516SKenneth E. Jansenc
46*59599516SKenneth E. Jansenc.... loop through the integration points
47*59599516SKenneth E. Jansenc
48*59599516SKenneth E. Jansen      rminv = zero
49*59599516SKenneth E. Jansen      rmass = zero
50*59599516SKenneth E. Jansen      qrl   = zero
51*59599516SKenneth E. Jansen
52*59599516SKenneth E. Jansen      do intp = 1, ngauss
53*59599516SKenneth E. Jansen
54*59599516SKenneth E. Jansen         call getshp(shp, shgl, sgn, shape, shdrv)
55*59599516SKenneth E. Jansen
56*59599516SKenneth E. Jansen         qdi = zero
57*59599516SKenneth E. Jansenc
58*59599516SKenneth E. Jansenc.... calculate the integration variables
59*59599516SKenneth E. Jansenc
60*59599516SKenneth E. Jansenc
61*59599516SKenneth E. Jansen         call e3qvar   (yl,           shdrv,
62*59599516SKenneth E. Jansen     &                  xl,           g1yi,
63*59599516SKenneth E. Jansen     &                  g2yi,         g3yi,         shg,
64*59599516SKenneth E. Jansen     &                  dxidx,        WdetJ )
65*59599516SKenneth E. Jansen
66*59599516SKenneth E. Jansen         call getdiff(dwl,  yl, shape, xmudmi, xl,rmu, tmp)
67*59599516SKenneth E. Jansenc
68*59599516SKenneth E. Jansenc.... diffusive flux in x1-direction
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansen         qdi(:,1) =  two * rmu *  g1yi(:,2)
71*59599516SKenneth E. Jansen         qdi(:,4) =        rmu * (g1yi(:,3) + g2yi(:,2))
72*59599516SKenneth E. Jansen         qdi(:,7) =        rmu * (g1yi(:,4) + g3yi(:,2))
73*59599516SKenneth E. Jansenc
74*59599516SKenneth E. Jansenc.... diffusive flux in x2-direction
75*59599516SKenneth E. Jansenc
76*59599516SKenneth E. Jansen         qdi(:,2) =        rmu * (g1yi(:,3) + g2yi(:,2))
77*59599516SKenneth E. Jansen         qdi(:,5) =  two * rmu *  g2yi(:,3)
78*59599516SKenneth E. Jansen         qdi(:,8) =        rmu * (g2yi(:,4) + g3yi(:,3))
79*59599516SKenneth E. Jansenc
80*59599516SKenneth E. Jansenc.... diffusive flux in x3-direction
81*59599516SKenneth E. Jansenc
82*59599516SKenneth E. Jansen         qdi(:,3) =        rmu * (g1yi(:,4) + g3yi(:,2))
83*59599516SKenneth E. Jansen         qdi(:,6)=        rmu * (g2yi(:,4) + g3yi(:,3))
84*59599516SKenneth E. Jansen         qdi(:,9)=  two * rmu *  g3yi(:,4)
85*59599516SKenneth E. Jansenc
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansenc.... assemble contribution of qdi to qrl,i.e., contribution to
88*59599516SKenneth E. Jansenc     each element shape function
89*59599516SKenneth E. Jansenc
90*59599516SKenneth E. Jansen         tmp = Qwt(lcsyst,intp)
91*59599516SKenneth E. Jansen         if (lcsyst .eq. 1) then
92*59599516SKenneth E. Jansen            tmp = tmp*(three/four)
93*59599516SKenneth E. Jansen         endif
94*59599516SKenneth E. Jansenc
95*59599516SKenneth E. Jansenc reconsider this when hierarchic wedges come into code WDGCHECK
96*59599516SKenneth E. Jansenc
97*59599516SKenneth E. Jansen
98*59599516SKenneth E. Jansen         do i=1,nshl
99*59599516SKenneth E. Jansen            qrl(:,i,1 ) = qrl(:,i,1 )+ shape(:,i)*tmp*qdi(:,1 )
100*59599516SKenneth E. Jansen            qrl(:,i,2 ) = qrl(:,i,2 )+ shape(:,i)*tmp*qdi(:,2 )
101*59599516SKenneth E. Jansen            qrl(:,i,3 ) = qrl(:,i,3 )+ shape(:,i)*tmp*qdi(:,3 )
102*59599516SKenneth E. Jansen
103*59599516SKenneth E. Jansen            qrl(:,i,4 ) = qrl(:,i,4 )+ shape(:,i)*tmp*qdi(:,4 )
104*59599516SKenneth E. Jansen            qrl(:,i,5 ) = qrl(:,i,5 )+ shape(:,i)*tmp*qdi(:,5 )
105*59599516SKenneth E. Jansen            qrl(:,i,6 ) = qrl(:,i,6 )+ shape(:,i)*tmp*qdi(:,6 )
106*59599516SKenneth E. Jansen
107*59599516SKenneth E. Jansen            qrl(:,i,7 ) = qrl(:,i,7 )+ shape(:,i)*tmp*qdi(:,7 )
108*59599516SKenneth E. Jansen            qrl(:,i,8 ) = qrl(:,i,8 )+ shape(:,i)*tmp*qdi(:,8 )
109*59599516SKenneth E. Jansen            qrl(:,i,9 ) = qrl(:,i,9 )+ shape(:,i)*tmp*qdi(:,9 )
110*59599516SKenneth E. Jansen         enddo
111*59599516SKenneth E. Jansenc
112*59599516SKenneth E. Jansenc.... add contribution to local mass matrix
113*59599516SKenneth E. Jansenc
114*59599516SKenneth E. Jansen
115*59599516SKenneth E. Jansen         if (have_local_mass .eq. 0) then
116*59599516SKenneth E. Jansen            do i=1,nshl
117*59599516SKenneth E. Jansen               do j=1,nshl
118*59599516SKenneth E. Jansen                  rmass(:,i,j) = rmass(:,i,j)+shape(:,i)*shape(:,j)*tmp
119*59599516SKenneth E. Jansen              enddo
120*59599516SKenneth E. Jansen           enddo
121*59599516SKenneth E. Jansen        endif
122*59599516SKenneth E. Jansenc
123*59599516SKenneth E. Jansenc.... end of the loop over integration points
124*59599516SKenneth E. Jansenc
125*59599516SKenneth E. Jansen      enddo
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansenc
128*59599516SKenneth E. Jansenc.... find the inverse of the local mass matrix for each element
129*59599516SKenneth E. Jansen
130*59599516SKenneth E. Jansen
131*59599516SKenneth E. Jansen         if (have_local_mass .eq. 0) then
132*59599516SKenneth E. Jansen            allocate (lmassinv(iblock)%p(npro,nshl,nshl))
133*59599516SKenneth E. Jansen
134*59599516SKenneth E. Jansen            do iel=1,npro
135*59599516SKenneth E. Jansen               do i=1,nshl      ! form the identy matrix
136*59599516SKenneth E. Jansen                  do j=1,nshl
137*59599516SKenneth E. Jansen                     lmassinv(iblock)%p(iel,i,j) = 0.0
138*59599516SKenneth E. Jansen                  enddo
139*59599516SKenneth E. Jansen                  lmassinv(iblock)%p(iel,i,i)=1.0
140*59599516SKenneth E. Jansen               enddo
141*59599516SKenneth E. Jansenc
142*59599516SKenneth E. Jansenc.... LU factor the mass matrix
143*59599516SKenneth E. Jansenc
144*59599516SKenneth E. Jansen               call ludcmp(rmass(iel,:,:),nshl,nshl,indx,d)
145*59599516SKenneth E. Jansenc
146*59599516SKenneth E. Jansenc.... back substitute with the identy matrix to find the
147*59599516SKenneth E. Jansenc     matrix inverse
148*59599516SKenneth E. Jansenc
149*59599516SKenneth E. Jansen               do j=1,nshl
150*59599516SKenneth E. Jansen                  call lubksb(rmass(iel,:,:),nshl,nshl,indx,
151*59599516SKenneth E. Jansen     &                        lmassinv(iblock)%p(iel,:,j))
152*59599516SKenneth E. Jansen               enddo
153*59599516SKenneth E. Jansen            enddo
154*59599516SKenneth E. Jansen            rminv(:,:,:) = lmassinv(iblock)%p(:,:,:)
155*59599516SKenneth E. Jansen         else
156*59599516SKenneth E. Jansen            rminv(:,:,:) = lmassinv(iblock)%p(:,:,:)
157*59599516SKenneth E. Jansen         endif
158*59599516SKenneth E. Jansenc
159*59599516SKenneth E. Jansenc.... find the modal coefficients of ql by multiplying by the inverse of
160*59599516SKenneth E. Jansenc     the local mass matrix
161*59599516SKenneth E. Jansenc
162*59599516SKenneth E. Jansen      do iel=1,npro
163*59599516SKenneth E. Jansen        do j=1,9
164*59599516SKenneth E. Jansenc         do j=1, 3*nsd
165*59599516SKenneth E. Jansen            ql(iel,:,j) = matmul( rminv(iel,:,:),qrl(iel,:,j) )
166*59599516SKenneth E. Jansen         enddo
167*59599516SKenneth E. Jansen      enddo
168*59599516SKenneth E. Jansenc
169*59599516SKenneth E. Jansenc.... return
170*59599516SKenneth E. Jansenc
171*59599516SKenneth E. Jansen      return
172*59599516SKenneth E. Jansen      end
173*59599516SKenneth E. Jansen
174*59599516SKenneth E. Jansen
175*59599516SKenneth E. Jansen
176*59599516SKenneth E. Jansen
177*59599516SKenneth E. Jansen      subroutine e3qlSclr (yl,      dwl,     shp,     shgl,
178*59599516SKenneth E. Jansen     &                     xl,      ql,      sgn )
179*59599516SKenneth E. Jansenc
180*59599516SKenneth E. Jansenc----------------------------------------------------------------------
181*59599516SKenneth E. Jansenc
182*59599516SKenneth E. Jansenc This routine computes the local diffusive flux vector using a
183*59599516SKenneth E. Jansenc local projection algorithm:
184*59599516SKenneth E. Jansenc     diffus * phi,i
185*59599516SKenneth E. Jansenc
186*59599516SKenneth E. Jansenc----------------------------------------------------------------------
187*59599516SKenneth E. Jansenc
188*59599516SKenneth E. Jansen      use local_mass
189*59599516SKenneth E. Jansen      include "common.h"
190*59599516SKenneth E. Jansenc
191*59599516SKenneth E. Jansen      dimension yl(npro,nshl,ndof),        dwl(npro,nshl),
192*59599516SKenneth E. Jansen     &          shp(nshl,ngauss),          shgl(nsd,nshl,ngauss),
193*59599516SKenneth E. Jansen     &          xl(npro,nenl,nsd),         sgn(npro,nshl),
194*59599516SKenneth E. Jansen     &          ql(npro,nshl,nsd)
195*59599516SKenneth E. Jansenc
196*59599516SKenneth E. Jansenc local arrays
197*59599516SKenneth E. Jansenc
198*59599516SKenneth E. Jansen      dimension dxidx(npro,nsd,nsd),       WdetJ(npro),
199*59599516SKenneth E. Jansen     &          diffus(npro),
200*59599516SKenneth E. Jansen     &          rminv(npro,nshl,nshl),
201*59599516SKenneth E. Jansen     &          qrl(npro,nshl,nsd)
202*59599516SKenneth E. Jansenc
203*59599516SKenneth E. Jansen      dimension qdi(npro,nsd),    shape(npro,nshl),
204*59599516SKenneth E. Jansen     &          shdrv(npro,nsd,nshl),      indx(nshl),
205*59599516SKenneth E. Jansen     &          rmass(npro,nshl,nshl),     gradT(npro,nsd),
206*59599516SKenneth E. Jansen     &          eviscv(npro)
207*59599516SKenneth E. Jansen
208*59599516SKenneth E. Jansenc
209*59599516SKenneth E. Jansenc.... loop through the integration points
210*59599516SKenneth E. Jansenc
211*59599516SKenneth E. Jansen      rminv = zero
212*59599516SKenneth E. Jansen      rmass = zero
213*59599516SKenneth E. Jansen      qrl   = zero
214*59599516SKenneth E. Jansen
215*59599516SKenneth E. Jansen      do intp = 1, ngauss
216*59599516SKenneth E. Jansen
217*59599516SKenneth E. Jansen         call getshp(shp, shgl, sgn, shape, shdrv)
218*59599516SKenneth E. Jansen
219*59599516SKenneth E. Jansen         qdi = zero
220*59599516SKenneth E. Jansenc
221*59599516SKenneth E. Jansenc.... calculate the integration variables
222*59599516SKenneth E. Jansenc
223*59599516SKenneth E. Jansenc
224*59599516SKenneth E. Jansen         call e3qvarSclr  (yl,           shdrv,        xl,
225*59599516SKenneth E. Jansen     &                     gradT,        dxidx,        WdetJ )
226*59599516SKenneth E. Jansenc
227*59599516SKenneth E. Jansenc....  call function to sort out diffusivity (at end of this file)
228*59599516SKenneth E. Jansenc
229*59599516SKenneth E. Jansen         call getdiffsclr(dwl,shape,yl, diffus)
230*59599516SKenneth E. Jansenc
231*59599516SKenneth E. Jansenc.... diffusive flux in x1-direction
232*59599516SKenneth E. Jansenc
233*59599516SKenneth E. Jansen         qdi(:,1) =  diffus * gradT(:,1)
234*59599516SKenneth E. Jansen         qdi(:,2) =  diffus * gradT(:,2)
235*59599516SKenneth E. Jansen         qdi(:,3) =  diffus * gradT(:,3)
236*59599516SKenneth E. Jansen
237*59599516SKenneth E. Jansenc
238*59599516SKenneth E. Jansenc.... assemble contribution of qdi to qrl,i.e., contribution to
239*59599516SKenneth E. Jansenc     each element shape function
240*59599516SKenneth E. Jansenc
241*59599516SKenneth E. Jansen         tmp = Qwt(lcsyst,intp)
242*59599516SKenneth E. Jansen         if (lcsyst .eq. 1) then
243*59599516SKenneth E. Jansen            tmp = tmp*(three/four)
244*59599516SKenneth E. Jansen         endif
245*59599516SKenneth E. Jansen
246*59599516SKenneth E. Jansen         do i=1,nshl
247*59599516SKenneth E. Jansen            qrl(:,i,1 ) = qrl(:,i,1 )+ shape(:,i)*tmp*qdi(:,1 )
248*59599516SKenneth E. Jansen            qrl(:,i,2 ) = qrl(:,i,2 )+ shape(:,i)*tmp*qdi(:,2 )
249*59599516SKenneth E. Jansen            qrl(:,i,3 ) = qrl(:,i,3 )+ shape(:,i)*tmp*qdi(:,3 )
250*59599516SKenneth E. Jansen         enddo
251*59599516SKenneth E. Jansenc
252*59599516SKenneth E. Jansenc.... add contribution to local mass matrix
253*59599516SKenneth E. Jansenc
254*59599516SKenneth E. Jansen         if (have_local_mass .eq. 0) then
255*59599516SKenneth E. Jansen            do i=1,nshl
256*59599516SKenneth E. Jansen               do j=1,nshl
257*59599516SKenneth E. Jansen                  rmass(:,i,j)=rmass(:,i,j)+shape(:,i)*shape(:,j)*tmp
258*59599516SKenneth E. Jansen               enddo
259*59599516SKenneth E. Jansen            enddo
260*59599516SKenneth E. Jansen         endif
261*59599516SKenneth E. Jansen
262*59599516SKenneth E. Jansenc.... end of the loop over integration points
263*59599516SKenneth E. Jansenc
264*59599516SKenneth E. Jansen      enddo
265*59599516SKenneth E. Jansen
266*59599516SKenneth E. Jansenc
267*59599516SKenneth E. Jansenc.... find the inverse of the local mass matrix for each element
268*59599516SKenneth E. Jansenc
269*59599516SKenneth E. Jansen       qrl   = qrl/6.d0
270*59599516SKenneth E. Jansenc
271*59599516SKenneth E. Jansenc.... Assuming that lmassinv was already computed for flow equations
272*59599516SKenneth E. Jansenc
273*59599516SKenneth E. Jansen       rmass = rmass/6.0
274*59599516SKenneth E. Jansenc
275*59599516SKenneth E. Jansenc.... for cubics, it cannot be precomputed, so compute and
276*59599516SKenneth E. Jansenc     save it the first time it is needed
277*59599516SKenneth E. Jansenc
278*59599516SKenneth E. Jansen         if (have_local_mass .eq. 0) then
279*59599516SKenneth E. Jansen            allocate (lmassinv(iblock)%p(npro,nshl,nshl))
280*59599516SKenneth E. Jansen
281*59599516SKenneth E. Jansen            do iel=1,npro
282*59599516SKenneth E. Jansen               do i=1,nshl      ! form the identy matrix
283*59599516SKenneth E. Jansen                  do j=1,nshl
284*59599516SKenneth E. Jansen                     lmassinv(iblock)%p(iel,i,j) = 0.0
285*59599516SKenneth E. Jansen                  enddo
286*59599516SKenneth E. Jansen                  lmassinv(iblock)%p(iel,i,i)=1.0
287*59599516SKenneth E. Jansen               enddo
288*59599516SKenneth E. Jansenc
289*59599516SKenneth E. Jansenc.... LU factor the mass matrix
290*59599516SKenneth E. Jansenc
291*59599516SKenneth E. Jansen               call ludcmp(rmass(iel,:,:),nshl,nshl,indx,d)
292*59599516SKenneth E. Jansenc
293*59599516SKenneth E. Jansenc.... back substitute with the identy matrix to find the
294*59599516SKenneth E. Jansenc     matrix inverse
295*59599516SKenneth E. Jansenc
296*59599516SKenneth E. Jansen               do j=1,nshl
297*59599516SKenneth E. Jansen                  call lubksb(rmass(iel,:,:),nshl,nshl,indx,
298*59599516SKenneth E. Jansen     &                        lmassinv(iblock)%p(iel,:,j))
299*59599516SKenneth E. Jansen               enddo
300*59599516SKenneth E. Jansen            enddo
301*59599516SKenneth E. Jansen            rminv(:,:,:) = lmassinv(iblock)%p(:,:,:)
302*59599516SKenneth E. Jansen         else
303*59599516SKenneth E. Jansen            rminv(:,:,:) = lmassinv(iblock)%p(:,:,:)
304*59599516SKenneth E. Jansen         endif
305*59599516SKenneth E. Jansenc
306*59599516SKenneth E. Jansenc.... find the modal coefficients of ql by multiplying by the inverse of
307*59599516SKenneth E. Jansenc     the local mass matrix
308*59599516SKenneth E. Jansenc
309*59599516SKenneth E. Jansen      do iel=1,npro
310*59599516SKenneth E. Jansen         do j=1,nsd
311*59599516SKenneth E. Jansen            ql(iel,:,j) = matmul( rminv(iel,:,:),qrl(iel,:,j) )
312*59599516SKenneth E. Jansen         enddo
313*59599516SKenneth E. Jansen      enddo
314*59599516SKenneth E. Jansenc
315*59599516SKenneth E. Jansenc.... return
316*59599516SKenneth E. Jansenc
317*59599516SKenneth E. Jansen      return
318*59599516SKenneth E. Jansen      end
319