xref: /phasta/phSolver/incompressible/e3q.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine e3q (yl,      dwl,     shp,     shgl,
2*59599516SKenneth E. Jansen     &                  xl,      ql,      rmassl,
3*59599516SKenneth E. Jansen     &                  xmudmi,  sgn )
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansenc----------------------------------------------------------------------
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc This routine computes the element contribution to the
8*59599516SKenneth E. Jansenc diffusive flux vector and the lumped mass matrix.
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,nshl,nsd)        : nodal coordinates at current step
15*59599516SKenneth E. Jansenc
16*59599516SKenneth E. Jansenc output:
17*59599516SKenneth E. Jansenc  ql     (npro,nshl,idflx) : element RHS diffusion residual
18*59599516SKenneth E. Jansenc  rmassl     (npro,nshl)        : element lumped mass matrix
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 yl(npro,nshl,ndof),     dwl(npro,nenl),
25*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),      shgl(nsd,nshl,ngauss),
26*59599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),
27*59599516SKenneth E. Jansen     &            ql(npro,nshl,idflx),  rmassl(npro,nshl),
28*59599516SKenneth E. Jansen     &            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     &            rmu(npro)
36*59599516SKenneth E. Jansenc
37*59599516SKenneth E. Jansen        dimension qdi(npro,idflx),alph1(npro),alph2(npro)
38*59599516SKenneth E. Jansenc
39*59599516SKenneth E. Jansen        dimension sgn(npro,nshl),          shape(npro,nshl),
40*59599516SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),    shpsum(npro)
41*59599516SKenneth E. Jansen
42*59599516SKenneth E. Jansen        real*8 tmp(npro)
43*59599516SKenneth E. Jansenc
44*59599516SKenneth E. Jansenc.... for surface tension
45*59599516SKenneth E. Jansenc
46*59599516SKenneth E. Jansen        dimension g1yti(npro),          g2yti(npro),
47*59599516SKenneth E. Jansen     &            g3yti(npro)
48*59599516SKenneth E. Jansen        integer idflow
49*59599516SKenneth E. Jansenc
50*59599516SKenneth E. Jansenc.... loop through the integration points
51*59599516SKenneth E. Jansenc
52*59599516SKenneth E. Jansen
53*59599516SKenneth E. Jansen
54*59599516SKenneth E. Jansen        alph1 = 0.d0
55*59599516SKenneth E. Jansen        alph2 = 0.d0
56*59599516SKenneth E. Jansen
57*59599516SKenneth E. Jansen        do intp = 1, ngauss
58*59599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
59*59599516SKenneth E. Jansenc
60*59599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
61*59599516SKenneth E. Jansen     &              shape,        shdrv)
62*59599516SKenneth E. Jansen
63*59599516SKenneth E. Jansenc
64*59599516SKenneth E. Jansenc.... initialize
65*59599516SKenneth E. Jansenc
66*59599516SKenneth E. Jansen        qdi = zero
67*59599516SKenneth E. Jansenc
68*59599516SKenneth E. Jansenc
69*59599516SKenneth E. Jansenc.... calculate the integration variables necessary for the
70*59599516SKenneth E. Jansenc     formation of q
71*59599516SKenneth E. Jansenc
72*59599516SKenneth E. Jansen
73*59599516SKenneth E. Jansen        call e3qvar   (yl,        shdrv,
74*59599516SKenneth E. Jansen     &                 xl,           g1yi,
75*59599516SKenneth E. Jansen     &                 g2yi,      g3yi,         shg,
76*59599516SKenneth E. Jansen     &                 dxidx,     WdetJ )
77*59599516SKenneth E. Jansenc
78*59599516SKenneth E. Jansen        idflow = 9   ! we ALWAYS save space for tau_{ij} in q_i
79*59599516SKenneth E. Jansen                     ! even if idiff is not greater than 1
80*59599516SKenneth E. Jansen
81*59599516SKenneth E. Jansen        if(idiff >= 1) then   !so taking care of all the idiff=1,3
82*59599516SKenneth E. Jansenc
83*59599516SKenneth E. Jansenc.... compute diffusive fluxes
84*59599516SKenneth E. Jansenc
85*59599516SKenneth E. Jansenc.... compute the viscosity
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansen        call getdiff(dwl, yl, shape, xmudmi, xl, rmu, tmp)
88*59599516SKenneth E. Jansenc
89*59599516SKenneth E. Jansenc.... diffusive flux in x1-direction
90*59599516SKenneth E. Jansenc
91*59599516SKenneth E. Jansen        qdi(:,1) =  two * rmu *  g1yi(:,2)
92*59599516SKenneth E. Jansen        qdi(:,4) =        rmu * (g1yi(:,3) + g2yi(:,2))
93*59599516SKenneth E. Jansen        qdi(:,7) =        rmu * (g1yi(:,4) + g3yi(:,2))
94*59599516SKenneth E. Jansenc
95*59599516SKenneth E. Jansenc.... diffusive flux in x2-direction
96*59599516SKenneth E. Jansenc
97*59599516SKenneth E. Jansen        qdi(:,2) =        rmu * (g1yi(:,3) + g2yi(:,2))
98*59599516SKenneth E. Jansen        qdi(:,5) =  two * rmu *  g2yi(:,3)
99*59599516SKenneth E. Jansen        qdi(:,8) =        rmu * (g2yi(:,4) + g3yi(:,3))
100*59599516SKenneth E. Jansenc
101*59599516SKenneth E. Jansenc.... diffusive flux in x3-direction
102*59599516SKenneth E. Jansenc
103*59599516SKenneth E. Jansen        qdi(:,3) =        rmu * (g1yi(:,4) + g3yi(:,2))
104*59599516SKenneth E. Jansen        qdi(:,6)=         rmu * (g2yi(:,4) + g3yi(:,3))
105*59599516SKenneth E. Jansen        qdi(:,9)=  two * rmu *  g3yi(:,4)
106*59599516SKenneth E. Jansenc
107*59599516SKenneth E. Jansenc
108*59599516SKenneth E. Jansenc.... assemble contribution of qdi to ql,i.e., contribution to
109*59599516SKenneth E. Jansenc     each element node
110*59599516SKenneth E. Jansenc
111*59599516SKenneth E. Jansen        do i=1,nshl
112*59599516SKenneth E. Jansen           ql(:,i,1 ) = ql(:,i,1 )+ shape(:,i)*WdetJ*qdi(:,1 )
113*59599516SKenneth E. Jansen           ql(:,i,2 ) = ql(:,i,2 )+ shape(:,i)*WdetJ*qdi(:,2 )
114*59599516SKenneth E. Jansen           ql(:,i,3 ) = ql(:,i,3 )+ shape(:,i)*WdetJ*qdi(:,3 )
115*59599516SKenneth E. Jansen
116*59599516SKenneth E. Jansen           ql(:,i,4 ) = ql(:,i,4 )+ shape(:,i)*WdetJ*qdi(:,4 )
117*59599516SKenneth E. Jansen           ql(:,i,5 ) = ql(:,i,5 )+ shape(:,i)*WdetJ*qdi(:,5 )
118*59599516SKenneth E. Jansen           ql(:,i,6 ) = ql(:,i,6 )+ shape(:,i)*WdetJ*qdi(:,6 )
119*59599516SKenneth E. Jansen
120*59599516SKenneth E. Jansen           ql(:,i,7 ) = ql(:,i,7 )+ shape(:,i)*WdetJ*qdi(:,7 )
121*59599516SKenneth E. Jansen           ql(:,i,8 ) = ql(:,i,8 )+ shape(:,i)*WdetJ*qdi(:,8 )
122*59599516SKenneth E. Jansen           ql(:,i,9 ) = ql(:,i,9 )+ shape(:,i)*WdetJ*qdi(:,9 )
123*59599516SKenneth E. Jansen
124*59599516SKenneth E. Jansen        enddo
125*59599516SKenneth E. Jansenc
126*59599516SKenneth E. Jansenc.... compute and assemble the element contribution to the lumped
127*59599516SKenneth E. Jansenc     mass matrix
128*59599516SKenneth E. Jansenc
129*59599516SKenneth E. Jansenc
130*59599516SKenneth E. Jansenc.... row sum technique
131*59599516SKenneth E. Jansenc
132*59599516SKenneth E. Jansen        if ( idiff == 1 ) then
133*59599516SKenneth E. Jansen           do i=1,nshl
134*59599516SKenneth E. Jansen              rmassl(:,i) = rmassl(:,i) + shape(:,i)*WdetJ
135*59599516SKenneth E. Jansen           enddo
136*59599516SKenneth E. Jansen        endif
137*59599516SKenneth E. Jansenc
138*59599516SKenneth E. Jansenc.... "special lumping technique" (Hughes p. 445)
139*59599516SKenneth E. Jansenc
140*59599516SKenneth E. Jansen        if ( idiff == 3 ) then
141*59599516SKenneth E. Jansen           shpsum = zero
142*59599516SKenneth E. Jansen           do i=1,nshl
143*59599516SKenneth E. Jansen              shpsum = shpsum + shape(:,i)*shape(:,i)
144*59599516SKenneth E. Jansen              rmassl(:,i)=rmassl(:,i)+shape(:,i)*shape(:,i)*WdetJ
145*59599516SKenneth E. Jansen           enddo
146*59599516SKenneth E. Jansen           alph1 = alph1+WdetJ
147*59599516SKenneth E. Jansen           alph2 = alph2+shpsum*WdetJ
148*59599516SKenneth E. Jansen        endif
149*59599516SKenneth E. Jansen      endif                     ! end of idiff=1 .or. 3
150*59599516SKenneth E. Jansenc
151*59599516SKenneth E. Jansen      if(isurf .eq. 1) then
152*59599516SKenneth E. Jansenc
153*59599516SKenneth E. Jansenc.... initialize
154*59599516SKenneth E. Jansenc
155*59599516SKenneth E. Jansen        g1yti   = zero
156*59599516SKenneth E. Jansen        g2yti   = zero
157*59599516SKenneth E. Jansen        g3yti   = zero
158*59599516SKenneth E. Jansenc
159*59599516SKenneth E. Jansenc.... calculate the integration variables necessary for the
160*59599516SKenneth E. Jansenc     formation of q
161*59599516SKenneth E. Jansenc
162*59599516SKenneth E. Jansenc.... compute the global gradient of Yt-variables, assuming 6th entry as
163*59599516SKenneth E. Jansenc.... the phase indicator function
164*59599516SKenneth E. Jansenc
165*59599516SKenneth E. Jansenc  Yt_{,x_i}=SUM_{a=1}^nshl (N_{a,x_i}(int) Yta)
166*59599516SKenneth E. Jansenc
167*59599516SKenneth E. Jansen        do n = 1, nshl
168*59599516SKenneth E. Jansen          g1yti(:)  = g1yti(:)  + shg(:,n,1) * yl(:,n,6)
169*59599516SKenneth E. Jansen          g2yti(:)  = g2yti(:)  + shg(:,n,2) * yl(:,n,6)
170*59599516SKenneth E. Jansen          g3yti(:)  = g3yti(:)  + shg(:,n,3) * yl(:,n,6)
171*59599516SKenneth E. Jansen        enddo
172*59599516SKenneth E. Jansenc
173*59599516SKenneth E. Jansenc    computing N_{b}*N_{a,x_i)*yta*WdetJ
174*59599516SKenneth E. Jansenc
175*59599516SKenneth E. Jansen        do i=1,nshl
176*59599516SKenneth E. Jansen           ql(:,i,idflow+1)  = ql(:,i,idflow+1)
177*59599516SKenneth E. Jansen     &                       + shape(:,i)*WdetJ*g1yti
178*59599516SKenneth E. Jansen           ql(:,i,idflow+2)  = ql(:,i,idflow+2)
179*59599516SKenneth E. Jansen     &                       + shape(:,i)*WdetJ*g2yti
180*59599516SKenneth E. Jansen           ql(:,i,idflow+3)  = ql(:,i,idflow+3)
181*59599516SKenneth E. Jansen     &                       + shape(:,i)*WdetJ*g3yti
182*59599516SKenneth E. Jansen           rmassl(:,i) = rmassl(:,i) + shape(:,i)*WdetJ
183*59599516SKenneth E. Jansen        enddo
184*59599516SKenneth E. Jansen      endif  !end of the isurf
185*59599516SKenneth E. Jansenc
186*59599516SKenneth E. Jansenc.... end of the loop over integration points
187*59599516SKenneth E. Jansenc
188*59599516SKenneth E. Jansen      enddo
189*59599516SKenneth E. Jansenc
190*59599516SKenneth E. Jansenc.... normalize the mass matrix for idiff == 3
191*59599516SKenneth E. Jansenc
192*59599516SKenneth E. Jansen      if ( idiff == 3 ) then
193*59599516SKenneth E. Jansen         do i=1,nshl
194*59599516SKenneth E. Jansen            rmassl(:,i) = rmassl(:,i)*alph1/alph2
195*59599516SKenneth E. Jansen         enddo
196*59599516SKenneth E. Jansen      endif
197*59599516SKenneth E. Jansen
198*59599516SKenneth E. Jansen
199*59599516SKenneth E. Jansenc
200*59599516SKenneth E. Jansenc.... return
201*59599516SKenneth E. Jansenc
202*59599516SKenneth E. Jansen       return
203*59599516SKenneth E. Jansen       end
204*59599516SKenneth E. Jansen
205*59599516SKenneth E. Jansen
206*59599516SKenneth E. Jansen        subroutine e3qSclr (yl,      dwl,     shp,     shgl,
207*59599516SKenneth E. Jansen     &                      xl,      ql,      rmassl,
208*59599516SKenneth E. Jansen     &                      sgn )
209*59599516SKenneth E. Jansenc
210*59599516SKenneth E. Jansenc----------------------------------------------------------------------
211*59599516SKenneth E. Jansenc
212*59599516SKenneth E. Jansenc This routine computes the element contribution to the
213*59599516SKenneth E. Jansenc diffusive flux vector and the lumped mass matrix.
214*59599516SKenneth E. Jansenc
215*59599516SKenneth E. Jansenc----------------------------------------------------------------------
216*59599516SKenneth E. Jansenc
217*59599516SKenneth E. Jansen        include "common.h"
218*59599516SKenneth E. Jansenc
219*59599516SKenneth E. Jansen        dimension yl(npro,nshl,ndof),    dwl(npro,nshl),
220*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),      shgl(nsd,nshl,ngauss),
221*59599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),
222*59599516SKenneth E. Jansen     &            ql(npro,nshl,nsd),     rmassl(npro,nshl)
223*59599516SKenneth E. Jansenc
224*59599516SKenneth E. Jansenc local arrays
225*59599516SKenneth E. Jansenc
226*59599516SKenneth E. Jansen        dimension gradT(npro,nsd),
227*59599516SKenneth E. Jansen     &            dxidx(npro,nsd,nsd),       WdetJ(npro)
228*59599516SKenneth E. Jansenc
229*59599516SKenneth E. Jansen        dimension qdi(npro,nsd),alph1(npro),alph2(npro)
230*59599516SKenneth E. Jansenc
231*59599516SKenneth E. Jansen        dimension sgn(npro,nshl),          shape(npro,nshl),
232*59599516SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),    shpsum(npro)
233*59599516SKenneth E. Jansen
234*59599516SKenneth E. Jansen        real*8 diffus(npro)
235*59599516SKenneth E. Jansenc
236*59599516SKenneth E. Jansenc.... loop through the integration points
237*59599516SKenneth E. Jansenc
238*59599516SKenneth E. Jansen
239*59599516SKenneth E. Jansen
240*59599516SKenneth E. Jansen        alph1 = 0.d0
241*59599516SKenneth E. Jansen        alph2 = 0.d0
242*59599516SKenneth E. Jansen
243*59599516SKenneth E. Jansen        do intp = 1, ngauss
244*59599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
245*59599516SKenneth E. Jansenc
246*59599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
247*59599516SKenneth E. Jansen     &              shape,        shdrv)
248*59599516SKenneth E. Jansen
249*59599516SKenneth E. Jansenc
250*59599516SKenneth E. Jansenc.... initialize
251*59599516SKenneth E. Jansenc
252*59599516SKenneth E. Jansen        qdi = zero
253*59599516SKenneth E. Jansenc
254*59599516SKenneth E. Jansenc
255*59599516SKenneth E. Jansenc.... calculate the integration variables necessary for the
256*59599516SKenneth E. Jansenc     formation of q
257*59599516SKenneth E. Jansenc
258*59599516SKenneth E. Jansen        call e3qvarSclr   (yl,        shdrv,
259*59599516SKenneth E. Jansen     &                     xl,        gradT,
260*59599516SKenneth E. Jansen     &                     dxidx,     WdetJ )
261*59599516SKenneth E. Jansen
262*59599516SKenneth E. Jansenc
263*59599516SKenneth E. Jansenc.... compute diffusive flux vector at this integration point
264*59599516SKenneth E. Jansenc
265*59599516SKenneth E. Jansen        call getdiffsclr(shape, dwl, yl, diffus)
266*59599516SKenneth E. Jansen
267*59599516SKenneth E. Jansenc
268*59599516SKenneth E. Jansenc.... diffusive flux
269*59599516SKenneth E. Jansenc
270*59599516SKenneth E. Jansen        qdi(:,1) =  diffus * gradT(:,1)
271*59599516SKenneth E. Jansen        qdi(:,2) =  diffus * gradT(:,2)
272*59599516SKenneth E. Jansen        qdi(:,3) =  diffus * gradT(:,3)
273*59599516SKenneth E. Jansenc
274*59599516SKenneth E. Jansenc
275*59599516SKenneth E. Jansenc.... assemble contribution of qdi to ql,i.e., contribution to
276*59599516SKenneth E. Jansenc     each element node
277*59599516SKenneth E. Jansenc
278*59599516SKenneth E. Jansen        do i=1,nshl
279*59599516SKenneth E. Jansen           ql(:,i,1 ) = ql(:,i,1 )+ shape(:,i)*WdetJ*qdi(:,1 )
280*59599516SKenneth E. Jansen           ql(:,i,2 ) = ql(:,i,2 )+ shape(:,i)*WdetJ*qdi(:,2 )
281*59599516SKenneth E. Jansen           ql(:,i,3 ) = ql(:,i,3 )+ shape(:,i)*WdetJ*qdi(:,3 )
282*59599516SKenneth E. Jansen
283*59599516SKenneth E. Jansen        enddo
284*59599516SKenneth E. Jansenc
285*59599516SKenneth E. Jansenc.... compute and assemble the element contribution to the lumped
286*59599516SKenneth E. Jansenc     mass matrix
287*59599516SKenneth E. Jansenc
288*59599516SKenneth E. Jansenc
289*59599516SKenneth E. Jansenc.... row sum technique
290*59599516SKenneth E. Jansenc
291*59599516SKenneth E. Jansen        if ( idiff == 1 ) then
292*59599516SKenneth E. Jansen           do i=1,nshl
293*59599516SKenneth E. Jansen              rmassl(:,i) = rmassl(:,i) + shape(:,i)*WdetJ
294*59599516SKenneth E. Jansen           enddo
295*59599516SKenneth E. Jansen        endif
296*59599516SKenneth E. Jansenc
297*59599516SKenneth E. Jansenc.... "special lumping technique" (Hughes p. 445)
298*59599516SKenneth E. Jansenc
299*59599516SKenneth E. Jansen        if ( idiff == 3 ) then
300*59599516SKenneth E. Jansen           shpsum = zero
301*59599516SKenneth E. Jansen           do i=1,nshl
302*59599516SKenneth E. Jansen              shpsum = shpsum + shape(:,i)*shape(:,i)
303*59599516SKenneth E. Jansen              rmassl(:,i)=rmassl(:,i)+shape(:,i)*shape(:,i)*WdetJ
304*59599516SKenneth E. Jansen           enddo
305*59599516SKenneth E. Jansen           alph1 = alph1+WdetJ
306*59599516SKenneth E. Jansen           alph2 = alph2+shpsum*WdetJ
307*59599516SKenneth E. Jansen        endif
308*59599516SKenneth E. Jansenc
309*59599516SKenneth E. Jansenc.... end of the loop over integration points
310*59599516SKenneth E. Jansenc
311*59599516SKenneth E. Jansen      enddo
312*59599516SKenneth E. Jansenc
313*59599516SKenneth E. Jansenc.... normalize the mass matrix for idiff == 3
314*59599516SKenneth E. Jansenc
315*59599516SKenneth E. Jansen      if ( idiff == 3 ) then
316*59599516SKenneth E. Jansen         do i=1,nshl
317*59599516SKenneth E. Jansen            rmassl(:,i) = rmassl(:,i)*alph1/alph2
318*59599516SKenneth E. Jansen         enddo
319*59599516SKenneth E. Jansen      endif
320*59599516SKenneth E. Jansenc
321*59599516SKenneth E. Jansenc.... return
322*59599516SKenneth E. Jansenc
323*59599516SKenneth E. Jansen       return
324*59599516SKenneth E. Jansen       end
325*59599516SKenneth E. Jansen
326