xref: /phasta/phSolver/incompressible/e3.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine e3 (yl,      acl,     dwl,     shp,
2*59599516SKenneth E. Jansen     &                 shgl,    xl,      rl,      ql,
3*59599516SKenneth E. Jansen     &                 xKebe,   xGoC,    xmudmi,  sgn,
4*59599516SKenneth E. Jansen     &                 rerrl, rlsl)
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc----------------------------------------------------------------------
7*59599516SKenneth E. Jansenc
8*59599516SKenneth E. Jansenc     This routine calculates the residual and tangent matrix for the
9*59599516SKenneth E. Jansenc     UBar formulation of the incompressible Navier Stokes equations.
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansenc
12*59599516SKenneth E. Jansenc input:    e    a   1..5   when we think of U^e_a  and U is 5 variables
13*59599516SKenneth E. Jansenc  yl     (npro,nshl,ndof)      : Y variables (not U)
14*59599516SKenneth E. Jansenc  acl    (npro,nshl,ndof)      : Y acceleration (Y_{,t})
15*59599516SKenneth E. Jansenc  shp    (nen,ngauss)           : element shape-functions  N_a
16*59599516SKenneth E. Jansenc  shgl   (nsd,nen,ngauss)       : element local-shape-functions N_{a,xi}
17*59599516SKenneth E. Jansenc  wght   (ngauss)               : element weight (for quadrature)
18*59599516SKenneth E. Jansenc  xl     (npro,nenl,nsd)       : nodal coordinates at current step (x^e_a)
19*59599516SKenneth E. Jansenc  ql     (npro,nshl,nsd*nsd) : diffusive flux vector (don't worry)
20*59599516SKenneth E. Jansenc  rlsl   (npro,nshl,6)       : resolved Leonard stresses
21*59599516SKenneth E. Jansenc
22*59599516SKenneth E. Jansenc output:
23*59599516SKenneth E. Jansenc  rl     (npro,nshl,nflow)      : element RHS residual    (G^e_a)
24*59599516SKenneth E. Jansenc  rml    (npro,nshl,nflow)      : element modified residual  (G^e_a tilde)
25*59599516SKenneth E. Jansenc  xKebe  (npro,9,nshl,nshl)  : element LHS tangent mass matrix
26*59599516SKenneth E. Jansenc  xGoC   (npro,4,nshl,nshl)    : element LHS tangent mass matrix
27*59599516SKenneth E. Jansenc
28*59599516SKenneth E. Jansenc Note: This routine will calculate the element matrices for the
29*59599516SKenneth E. Jansenc        Hulbert's generalized alpha method integrator
30*59599516SKenneth E. Jansenc
31*59599516SKenneth E. Jansenc Mathematics done by:  Michel Mallet, Farzin Shakib (version 1)
32*59599516SKenneth E. Jansenc                       Farzin Shakib                (version 2)
33*59599516SKenneth E. Jansenc
34*59599516SKenneth E. Jansenc K. E. Jansen,   Winter 1999.   (advective form formulation)
35*59599516SKenneth E. Jansenc C. H. Whiting,  Winter 1999.   (advective form formulation)
36*59599516SKenneth E. Jansenc----------------------------------------------------------------------
37*59599516SKenneth E. Jansenc
38*59599516SKenneth E. Jansen        include "common.h"
39*59599516SKenneth E. Jansenc
40*59599516SKenneth E. Jansen        dimension yl(npro,nshl,ndof),
41*59599516SKenneth E. Jansen     &            acl(npro,nshl,ndof),
42*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),       shgl(nsd,nshl,ngauss),
43*59599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),      dwl(npro,nenl),
44*59599516SKenneth E. Jansen     &            rl(npro,nshl,nflow),     ql(npro,nshl,idflx)
45*59599516SKenneth E. Jansenc
46*59599516SKenneth E. Jansen
47*59599516SKenneth E. Jansen        dimension xKebe(npro,9,nshl,nshl), xGoC(npro,4,nshl,nshl)
48*59599516SKenneth E. Jansenc
49*59599516SKenneth E. Jansenc.... local declarations
50*59599516SKenneth E. Jansenc
51*59599516SKenneth E. Jansen        dimension g1yi(npro,ndof),        g2yi(npro,ndof),
52*59599516SKenneth E. Jansen     &            g3yi(npro,ndof),        shg(npro,nshl,nsd),
53*59599516SKenneth E. Jansen     &            aci(npro,3),            dxidx(npro,nsd,nsd),
54*59599516SKenneth E. Jansen     &            WdetJ(npro),            rho(npro),
55*59599516SKenneth E. Jansen     &            pres(npro),             u1(npro),
56*59599516SKenneth E. Jansen     &            u2(npro),               u3(npro),
57*59599516SKenneth E. Jansen     &            rLui(npro,nsd),         uBar(npro,nsd),
58*59599516SKenneth E. Jansen     &            xmudmi(npro,ngauss),     sgn(npro,nshl),
59*59599516SKenneth E. Jansen     &            shpfun(npro,nshl),      shdrv(npro,nsd,nshl),
60*59599516SKenneth E. Jansen     &            rmu(npro),              tauC(npro),
61*59599516SKenneth E. Jansen     &            tauM(npro),             tauBar(npro),
62*59599516SKenneth E. Jansen     &            src(npro,3)
63*59599516SKenneth E. Jansen
64*59599516SKenneth E. Jansen        dimension rlsl(npro,nshl,6),      rlsli(npro,6)
65*59599516SKenneth E. Jansen
66*59599516SKenneth E. Jansen        real*8    rerrl(npro,nshl,6)
67*59599516SKenneth E. Jansen        integer   aa
68*59599516SKenneth E. Jansen
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansenc
71*59599516SKenneth E. Jansenc.... local reconstruction of diffusive flux vector for quadratics
72*59599516SKenneth E. Jansenc     or greater but NOT for bflux since local mass was not mapped
73*59599516SKenneth E. Jansenc
74*59599516SKenneth E. Jansen        if ( idiff==2 .and. ires .eq. 1 ) then
75*59599516SKenneth E. Jansen           call e3ql (yl,        dwl,       shp,       shgl,
76*59599516SKenneth E. Jansen     &                xl,        ql,        xmudmi,
77*59599516SKenneth E. Jansen     &                sgn)
78*59599516SKenneth E. Jansen        endif
79*59599516SKenneth E. Jansenc
80*59599516SKenneth E. Jansenc.... loop through the integration points
81*59599516SKenneth E. Jansenc
82*59599516SKenneth E. Jansen
83*59599516SKenneth E. Jansen        do intp = 1, ngauss
84*59599516SKenneth E. Jansen
85*59599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point
88*59599516SKenneth E. Jansenc
89*59599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
90*59599516SKenneth E. Jansen     &              shpfun,       shdrv)
91*59599516SKenneth E. Jansenc
92*59599516SKenneth E. Jansenc.... get necessary fluid properties (including eddy viscosity)
93*59599516SKenneth E. Jansenc
94*59599516SKenneth E. Jansen        call getdiff(dwl,  yl,     shpfun,     xmudmi, xl,   rmu, rho)
95*59599516SKenneth E. Jansenc
96*59599516SKenneth E. Jansenc.... calculate the integration variables
97*59599516SKenneth E. Jansenc
98*59599516SKenneth E. Jansen        call e3ivar (yl,          acl,       shpfun,
99*59599516SKenneth E. Jansen     &               shdrv,       xl,
100*59599516SKenneth E. Jansen     &               aci,         g1yi,      g2yi,
101*59599516SKenneth E. Jansen     &               g3yi,        shg,       dxidx,
102*59599516SKenneth E. Jansen     &               WdetJ,       rho,       pres,
103*59599516SKenneth E. Jansen     &               u1,          u2,        u3,
104*59599516SKenneth E. Jansen     &               ql,          rLui,      src,
105*59599516SKenneth E. Jansen     &               rerrl,       rlsl,      rlsli,
106*59599516SKenneth E. Jansen     &               dwl)
107*59599516SKenneth E. Jansenc
108*59599516SKenneth E. Jansenc.... compute the stabilization terms
109*59599516SKenneth E. Jansenc
110*59599516SKenneth E. Jansen        call e3stab (rho,          u1,       u2,
111*59599516SKenneth E. Jansen     &               u3,           dxidx,    rLui,
112*59599516SKenneth E. Jansen     &               rmu,          tauC,     tauM,
113*59599516SKenneth E. Jansen     &               tauBar,       uBar )
114*59599516SKenneth E. Jansenc
115*59599516SKenneth E. Jansenc.... compute the residual contribution at this integration point
116*59599516SKenneth E. Jansenc
117*59599516SKenneth E. Jansen        call e3Res ( u1,        u2,         u3,
118*59599516SKenneth E. Jansen     &               uBar,      aci,        WdetJ,
119*59599516SKenneth E. Jansen     &               g1yi,      g2yi,       g3yi,
120*59599516SKenneth E. Jansen     &               rLui,      rmu,        rho,
121*59599516SKenneth E. Jansen     &               tauC,      tauM,       tauBar,
122*59599516SKenneth E. Jansen     &               shpfun,    shg,        src,
123*59599516SKenneth E. Jansen     &               rl,        pres,       acl,
124*59599516SKenneth E. Jansen     &               rlsli)
125*59599516SKenneth E. Jansenc
126*59599516SKenneth E. Jansenc.... compute the tangent matrix contribution
127*59599516SKenneth E. Jansenc
128*59599516SKenneth E. Jansen        if (lhs .eq. 1) then
129*59599516SKenneth E. Jansen           call e3LHS ( u1,        u2,         u3,
130*59599516SKenneth E. Jansen     &                  uBar,      WdetJ,      rho,
131*59599516SKenneth E. Jansen     &                  rLui,      rmu,
132*59599516SKenneth E. Jansen     &                  tauC,      tauM,       tauBar,
133*59599516SKenneth E. Jansen     &                  shpfun,    shg,        xKebe,
134*59599516SKenneth E. Jansen     &                  xGoC )
135*59599516SKenneth E. Jansen        endif
136*59599516SKenneth E. Jansen
137*59599516SKenneth E. Jansenc
138*59599516SKenneth E. Jansenc.... end of integration loop
139*59599516SKenneth E. Jansenc
140*59599516SKenneth E. Jansen      enddo
141*59599516SKenneth E. Jansen
142*59599516SKenneth E. Jansenc
143*59599516SKenneth E. Jansenc.... symmetrize C
144*59599516SKenneth E. Jansenc
145*59599516SKenneth E. Jansen      if (lhs .eq. 1) then
146*59599516SKenneth E. Jansen         do ib = 1, nshl
147*59599516SKenneth E. Jansen            do iaa = 1, ib-1
148*59599516SKenneth E. Jansen               xGoC(:,4,iaa,ib) = xGoC(:,4,ib,iaa)
149*59599516SKenneth E. Jansen            enddo
150*59599516SKenneth E. Jansen         enddo
151*59599516SKenneth E. Jansen      endif
152*59599516SKenneth E. Jansenc
153*59599516SKenneth E. Jansenc.... return
154*59599516SKenneth E. Jansenc
155*59599516SKenneth E. Jansen      return
156*59599516SKenneth E. Jansen      end
157*59599516SKenneth E. Jansen
158*59599516SKenneth E. Jansen
159*59599516SKenneth E. Jansenc^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
160*59599516SKenneth E. Jansenc###################################################################
161*59599516SKenneth E. Jansen
162*59599516SKenneth E. Jansen
163*59599516SKenneth E. Jansen      subroutine e3Sclr (yl,      acl,     shp,
164*59599516SKenneth E. Jansen     &                     shgl,    xl,      dwl,
165*59599516SKenneth E. Jansen     &                     rl,      ql,      xSebe,
166*59599516SKenneth E. Jansen     &                     sgn,     xmudmi)
167*59599516SKenneth E. Jansenc
168*59599516SKenneth E. Jansenc----------------------------------------------------------------------
169*59599516SKenneth E. Jansenc
170*59599516SKenneth E. Jansenc     This routine calculates the residual and tangent matrix for the
171*59599516SKenneth E. Jansenc     advection - diffusion equation for scalar.
172*59599516SKenneth E. Jansenc
173*59599516SKenneth E. Jansenc K. E. Jansen,   Winter 1999.   (advective form formulation)
174*59599516SKenneth E. Jansenc C. H. Whiting,  Winter 1999.   (advective form formulation)
175*59599516SKenneth E. Jansenc----------------------------------------------------------------------
176*59599516SKenneth E. Jansenc
177*59599516SKenneth E. Jansen      include "common.h"
178*59599516SKenneth E. Jansenc
179*59599516SKenneth E. Jansen      real*8    yl(npro,nshl,ndof),     acl(npro,nshl,ndof),
180*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),       shgl(nsd,nshl,ngauss),
181*59599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),      rl(npro,nshl),
182*59599516SKenneth E. Jansen     &            ql(npro,nshl,nsd),      xSebe(npro,nshl,nshl),
183*59599516SKenneth E. Jansen     &            dwl(npro,nshl)
184*59599516SKenneth E. Jansenc
185*59599516SKenneth E. Jansenc.... local declarations
186*59599516SKenneth E. Jansenc
187*59599516SKenneth E. Jansen      real*8    gradS(npro,nsd),        shg(npro,nshl,nsd),
188*59599516SKenneth E. Jansen     &            Sdot(npro),             Sclr(npro),
189*59599516SKenneth E. Jansen     &            dxidx(npro,nsd,nsd),    WdetJ(npro),
190*59599516SKenneth E. Jansen     &            u1(npro),     u2(npro), u3(npro),
191*59599516SKenneth E. Jansen     &            sgn(npro,nshl),         shpfun(npro,nshl),
192*59599516SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),   rLS(npro),
193*59599516SKenneth E. Jansen     &            tauS(npro),             diffus(npro),
194*59599516SKenneth E. Jansen     &            srcL(npro),             srcR(npro),
195*59599516SKenneth E. Jansen     &            gGradS(npro,nsd),       dcFct(npro),
196*59599516SKenneth E. Jansen     &            giju(npro,6)
197*59599516SKenneth E. Jansenc
198*59599516SKenneth E. Jansenc.... Source terms sometimes take the form (beta_i)*(phi,_i).  Since
199*59599516SKenneth E. Jansenc     the convective term has (u_i)*(phi,_i), it is useful to treat
200*59599516SKenneth E. Jansenc     beta_i as a "correction" to the velocity.  In calculating the
201*59599516SKenneth E. Jansenc     stabilization terms, the new "modified" velocity (u_i-beta_i) is
202*59599516SKenneth E. Jansenc     then used in place of the pure velocity for stabilization terms,
203*59599516SKenneth E. Jansenc     and the source term sneaks into the RHS and LHS.
204*59599516SKenneth E. Jansen      real*8    uMod(npro,nsd), srcRat(npro), xmudmi(npro,ngauss)
205*59599516SKenneth E. Jansenc
206*59599516SKenneth E. Jansen      integer   aa, b
207*59599516SKenneth E. Jansenc
208*59599516SKenneth E. Jansenc.... local reconstruction of diffusive flux vector
209*59599516SKenneth E. Jansenc
210*59599516SKenneth E. Jansen        if ( idiff==2 ) then
211*59599516SKenneth E. Jansen           call e3qlSclr (yl, dwl, shp, shgl, xl, ql, sgn)
212*59599516SKenneth E. Jansen        endif
213*59599516SKenneth E. Jansenc
214*59599516SKenneth E. Jansenc.... loop through the integration points
215*59599516SKenneth E. Jansenc
216*59599516SKenneth E. Jansen        do intp = 1, ngauss
217*59599516SKenneth E. Jansen
218*59599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
219*59599516SKenneth E. Jansenc
220*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point
221*59599516SKenneth E. Jansenc
222*59599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
223*59599516SKenneth E. Jansen     &              shpfun,        shdrv)
224*59599516SKenneth E. Jansenc
225*59599516SKenneth E. Jansenc.... get necessary fluid properties
226*59599516SKenneth E. Jansenc
227*59599516SKenneth E. Jansen        call getdiffsclr(shpfun,dwl,yl,diffus)
228*59599516SKenneth E. Jansenc
229*59599516SKenneth E. Jansenc.... calculate the integration variables
230*59599516SKenneth E. Jansenc
231*59599516SKenneth E. Jansen        call e3ivarSclr(yl,          acl,       shpfun,
232*59599516SKenneth E. Jansen     &                  shdrv,       xl,        xmudmi,
233*59599516SKenneth E. Jansen     &                  Sclr,        Sdot,      gradS,
234*59599516SKenneth E. Jansen     &                  shg,         dxidx,     WdetJ,
235*59599516SKenneth E. Jansen     &                  u1,          u2,        u3,
236*59599516SKenneth E. Jansen     &                  ql,          rLS,       SrcR,
237*59599516SKenneth E. Jansen     &                  SrcL,        uMod,      dwl,
238*59599516SKenneth E. Jansen     &                  diffus,      srcRat)
239*59599516SKenneth E. Jansen
240*59599516SKenneth E. Jansen
241*59599516SKenneth E. Jansenc
242*59599516SKenneth E. Jansenc.... compute the stabilization terms
243*59599516SKenneth E. Jansenc
244*59599516SKenneth E. Jansen        call e3StabSclr (uMod,    dxidx,   tauS,
245*59599516SKenneth E. Jansen     &                   diffus,  srcR,    giju,
246*59599516SKenneth E. Jansen     &                   srcRat)
247*59599516SKenneth E. Jansenc
248*59599516SKenneth E. Jansenc... computing the DC factor for the discontinuity capturing
249*59599516SKenneth E. Jansenc
250*59599516SKenneth E. Jansen        if (idcsclr(1) .ne. 0) then
251*59599516SKenneth E. Jansen           if ((idcsclr(2).eq.1 .and. isclr.eq.1) .or.
252*59599516SKenneth E. Jansen     &          (idcsclr(2).eq.2 .and. isclr.eq.2)) then ! scalar with dc
253*59599516SKenneth E. Jansenc
254*59599516SKenneth E. Jansen              call e3dcSclr ( gradS,    giju,     gGradS,
255*59599516SKenneth E. Jansen     &                        rLS,      tauS,     srcR,
256*59599516SKenneth E. Jansen     &                        dcFct)
257*59599516SKenneth E. Jansen           endif
258*59599516SKenneth E. Jansen        endif                   !end of idcsclr
259*59599516SKenneth E. Jansenc
260*59599516SKenneth E. Jansenc.... compute the residual contribution at this integration point
261*59599516SKenneth E. Jansenc
262*59599516SKenneth E. Jansen        call e3ResSclr ( uMod,      gGradS,
263*59599516SKenneth E. Jansen     &                   Sclr,      Sdot,       gradS,
264*59599516SKenneth E. Jansen     &                   WdetJ,     rLS,        tauS,
265*59599516SKenneth E. Jansen     &                   shpfun,    shg,        srcR,
266*59599516SKenneth E. Jansen     &                   diffus,
267*59599516SKenneth E. Jansen     &                   rl )
268*59599516SKenneth E. Jansenc
269*59599516SKenneth E. Jansenc.... compute the tangent matrix contribution
270*59599516SKenneth E. Jansenc
271*59599516SKenneth E. Jansen        if (lhs .eq. 1) then
272*59599516SKenneth E. Jansen           call e3LHSSclr ( uMod,      giju,       dcFct,
273*59599516SKenneth E. Jansen     &                      Sclr,      Sdot,       gradS,
274*59599516SKenneth E. Jansen     &                      WdetJ,     rLS,        tauS,
275*59599516SKenneth E. Jansen     &                      shpfun,    shg,        srcL,
276*59599516SKenneth E. Jansen     &                      diffus,
277*59599516SKenneth E. Jansen     &                      xSebe )
278*59599516SKenneth E. Jansen
279*59599516SKenneth E. Jansen        endif
280*59599516SKenneth E. Jansen
281*59599516SKenneth E. Jansenc
282*59599516SKenneth E. Jansenc.... end of integration loop
283*59599516SKenneth E. Jansenc
284*59599516SKenneth E. Jansen      enddo
285*59599516SKenneth E. Jansen
286*59599516SKenneth E. Jansenc
287*59599516SKenneth E. Jansenc.... return
288*59599516SKenneth E. Jansenc
289*59599516SKenneth E. Jansen      return
290*59599516SKenneth E. Jansen      end
291