xref: /phasta/phSolver/compressible/e3.f (revision f0b806cbcd176bc4c068569bf06644e2ab71f1d7)
159599516SKenneth E. Jansen        subroutine e3 (yl,      ycl,     acl,     shp,
259599516SKenneth E. Jansen     &                 shgl,    xl,      rl,      rml,    xmudmi,
359599516SKenneth E. Jansen     &                 BDiagl,  ql,      sgn,     rlsl,   EGmass,
459599516SKenneth E. Jansen     &                 rerrl,   ytargetl)
559599516SKenneth E. Jansenc
659599516SKenneth E. Jansenc----------------------------------------------------------------------
759599516SKenneth E. Jansenc
859599516SKenneth E. Jansenc This routine is the 3D element routine for the N-S equations.
959599516SKenneth E. Jansenc This routine calculates the RHS residual and if requested the
1059599516SKenneth E. Jansenc modified residual or the LHS consistent mass matrix or the block-
1159599516SKenneth E. Jansenc diagonal preconditioner.
1259599516SKenneth E. Jansenc
1359599516SKenneth E. Jansenc input:
1459599516SKenneth E. Jansenc  yl     (npro,nshl,nflow)     : Y variables  (DOES NOT CONTAIN SCALARS)
1559599516SKenneth E. Jansenc  ycl    (npro,nshl,ndof)      : Y variables at current step
1659599516SKenneth E. Jansenc  acl    (npro,nshl,ndof)      : Y acceleration (Y_{,t})
1759599516SKenneth E. Jansenc  shp    (nshl,ngauss)       : element shape-functions  N_a
1859599516SKenneth E. Jansenc  shgl   (nsd,nshl,ngauss)   : element local-shape-functions N_{a,xi}
1959599516SKenneth E. Jansenc  xl     (npro,nenl,nsd)       : nodal coordinates at current step (x^e_a)
2059599516SKenneth E. Jansenc  ql     (npro,nshl,(nflow-1)*nsd) : diffusive flux vector
2159599516SKenneth E. Jansenc  rlsl   (npro,nshl,6)       : resolved Leonard stresses
2259599516SKenneth E. Jansenc  sgn    (npro,nshl)         : shape function sign matrix
2359599516SKenneth E. Jansenc
2459599516SKenneth E. Jansenc output:
2559599516SKenneth E. Jansenc  rl     (npro,nshl,nflow)      : element RHS residual    (G^e_a)
2659599516SKenneth E. Jansenc  rml    (npro,nshl,nflow)      : element modified residual  (G^e_a tilde)
2759599516SKenneth E. Jansenc  EGmass (npro,nedof,nedof)    : element LHS tangent mass matrix (dG^e_a
2859599516SKenneth E. Jansenc                                                                  dY_b  )
2959599516SKenneth E. Jansenc  BDiagl (npro,nshl,nflow,nflow) : block-diagonal preconditioner
3059599516SKenneth E. Jansenc
3159599516SKenneth E. Jansenc
3259599516SKenneth E. Jansenc Note: This routine will calculate the element matrices for the
3359599516SKenneth E. Jansenc        Hulbert's generalized alpha method integrator
3459599516SKenneth E. Jansenc
3559599516SKenneth E. Jansenc Note: nedof = nflow*nshape is the total number of degrees of freedom
3659599516SKenneth E. Jansenc       at each element node.
3759599516SKenneth E. Jansenc
3859599516SKenneth E. Jansenc Mathematics done by:  Michel Mallet, Farzin Shakib (version 1)
3959599516SKenneth E. Jansenc                       Farzin Shakib                (version 2)
4059599516SKenneth E. Jansenc
4159599516SKenneth E. Jansenc
4259599516SKenneth E. Jansenc Zdenek Johan, Summer 1990.   (Modified from e2.f)
4359599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.   (Fortran 90)
4459599516SKenneth E. Jansenc Kenneth Jansen, Winter 1997. (Primitive Variables)
4559599516SKenneth E. Jansenc Chris Whiting, Winter 1998.  (LHS matrix formation)
4659599516SKenneth E. Jansenc----------------------------------------------------------------------
4759599516SKenneth E. Jansenc
4859599516SKenneth E. Jansen        include "common.h"
4959599516SKenneth E. Jansenc
5059599516SKenneth E. Jansen        dimension yl(npro,nshl,nflow),     ycl(npro,nshl,ndof),
5159599516SKenneth E. Jansen     &            acl(npro,nshl,ndof),     rmu(npro),
5259599516SKenneth E. Jansen     &            shp(nshl,ngauss),        rlm2mu(npro),
5359599516SKenneth E. Jansen     &            shgl(nsd,nshl,ngauss),   con(npro),
5459599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),       rlm(npro),
5559599516SKenneth E. Jansen     &            rl(npro,nshl,nflow),     ql(npro,nshl,idflx),
5659599516SKenneth E. Jansen     &            rml(npro,nshl,nflow),    xmudmi(npro,ngauss),
5759599516SKenneth E. Jansen     &            BDiagl(npro,nshl,nflow,nflow),
5859599516SKenneth E. Jansen     &            EGmass(npro,nedof,nedof),cv(npro),
5959599516SKenneth E. Jansen     &            ytargetl(npro,nshl,nflow)
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansen        dimension dui(npro,ndof),            aci(npro,ndof)
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansen        dimension g1yi(npro,nflow),          g2yi(npro,nflow),
6459599516SKenneth E. Jansen     &            g3yi(npro,nflow),          shg(npro,nshl,nsd),
6559599516SKenneth E. Jansen     &            divqi(npro,nflow),       tau(npro,5)
6659599516SKenneth E. Jansenc
6759599516SKenneth E. Jansen        dimension dxidx(npro,nsd,nsd),       WdetJ(npro)
6859599516SKenneth E. Jansenc
6959599516SKenneth E. Jansen        dimension rho(npro),                 pres(npro),
7059599516SKenneth E. Jansen     &            T(npro),                   ei(npro),
7159599516SKenneth E. Jansen     &            h(npro),                   alfap(npro),
7259599516SKenneth E. Jansen     &            betaT(npro),               DC(npro,ngauss),
7359599516SKenneth E. Jansen     &            cp(npro),                  rk(npro),
7459599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
7559599516SKenneth E. Jansen     &            u3(npro),                  A0DC(npro,4),
7659599516SKenneth E. Jansen     &            Sclr(npro),                dVdY(npro,15),
7759599516SKenneth E. Jansen     &            giju(npro,6),              rTLS(npro),
7859599516SKenneth E. Jansen     &            raLS(npro),                A0inv(npro,15)
7959599516SKenneth E. Jansenc
8059599516SKenneth E. Jansen        dimension A0(npro,nflow,nflow),      A1(npro,nflow,nflow),
8159599516SKenneth E. Jansen     &            A2(npro,nflow,nflow),      A3(npro,nflow,nflow)
8259599516SKenneth E. Jansenc
8359599516SKenneth E. Jansen        dimension rLyi(npro,nflow),          sgn(npro,nshl)
8459599516SKenneth E. Jansenc
8559599516SKenneth E. Jansen        dimension ri(npro,nflow*(nsd+1)),    rmi(npro,nflow*(nsd+1)),
8659599516SKenneth E. Jansen     &            shape(npro,nshl),          shdrv(npro,nsd,nshl),
8759599516SKenneth E. Jansen     &            stiff(npro,nsd*nflow,nsd*nflow),
8859599516SKenneth E. Jansen     &            PTau(npro,5,5),
8959599516SKenneth E. Jansen     &            sforce(npro,3),      compK(npro,10)
9059599516SKenneth E. Jansenc
9159599516SKenneth E. Jansen        dimension x(npro,3),              bcool(npro)
9259599516SKenneth E. Jansen
9359599516SKenneth E. Jansen        dimension rlsl(npro,nshl,6),      rlsli(npro,6)
9459599516SKenneth E. Jansen
9559599516SKenneth E. Jansen        real*8    rerrl(npro,nshl,6)
9659599516SKenneth E. Jansen        ttim(6) = ttim(6) - secs(0.0)
9759599516SKenneth E. Jansenc
9859599516SKenneth E. Jansenc.... local reconstruction of diffusive flux vector
9959599516SKenneth E. Jansenc (note: not currently included in mfg)
10059599516SKenneth E. Jansen        if (idiff==2 .and. (ires==3 .or. ires==1)) then
10159599516SKenneth E. Jansen           call e3ql (ycl, shp, shgl, xl, ql, xmudmi, sgn)
10259599516SKenneth E. Jansen        endif
10359599516SKenneth E. Jansenc
10459599516SKenneth E. Jansenc.... loop through the integration points
10559599516SKenneth E. Jansenc
10659599516SKenneth E. Jansen        do intp = 1, ngauss
10759599516SKenneth E. Jansenc
10859599516SKenneth E. Jansenc.... if Det. .eq. 0, do not include this point
10959599516SKenneth E. Jansenc
11059599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
11159599516SKenneth E. Jansenc
11259599516SKenneth E. Jansenc.... create a matrix of shape functions (and derivatives) for each
11359599516SKenneth E. Jansenc     element at this quadrature point. These arrays will contain
11459599516SKenneth E. Jansenc     the correct signs for the hierarchic basis
11559599516SKenneth E. Jansenc
11659599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
11759599516SKenneth E. Jansen     &              shape,        shdrv)
11859599516SKenneth E. Jansenc
11959599516SKenneth E. Jansenc.... initialize
12059599516SKenneth E. Jansenc
12159599516SKenneth E. Jansen        ri  = zero
12259599516SKenneth E. Jansen        rmi = zero
12359599516SKenneth E. Jansen        if (lhs .eq. 1) stiff = zero
12459599516SKenneth E. Jansenc
12559599516SKenneth E. Jansenc
12659599516SKenneth E. Jansenc.... calculate the integration variables
12759599516SKenneth E. Jansenc
12859599516SKenneth E. Jansen        ttim(8) = ttim(8) - secs(0.0)
12959599516SKenneth E. Jansen        call e3ivar (yl,              ycl,             acl,
13059599516SKenneth E. Jansen     &               Sclr,            shape,           shdrv,
13159599516SKenneth E. Jansen     &               xl,              dui,             aci,
13259599516SKenneth E. Jansen     &               g1yi,            g2yi,            g3yi,
13359599516SKenneth E. Jansen     &               shg,             dxidx,           WdetJ,
13459599516SKenneth E. Jansen     &               rho,             pres,            T,
13559599516SKenneth E. Jansen     &               ei,              h,               alfap,
13659599516SKenneth E. Jansen     &               betaT,           cp,              rk,
13759599516SKenneth E. Jansen     &               u1,              u2,              u3,
13859599516SKenneth E. Jansen     &               ql,              divqi,           sgn,
13959599516SKenneth E. Jansen     &               rLyi,  !passed as a work array
14059599516SKenneth E. Jansen     &               rmu,             rlm,             rlm2mu,
14159599516SKenneth E. Jansen     &               con,             rlsl,            rlsli,
14259599516SKenneth E. Jansen     &               xmudmi,          sforce,          cv)
14359599516SKenneth E. Jansen        ttim(8) = ttim(8) + secs(0.0)
14459599516SKenneth E. Jansen
14559599516SKenneth E. Jansenc
14659599516SKenneth E. Jansenc.... calculate the relevant matrices
14759599516SKenneth E. Jansenc
14859599516SKenneth E. Jansen        ttim(9) = ttim(9) - secs(0.0)
14959599516SKenneth E. Jansen        call e3mtrx (rho,             pres,           T,
15059599516SKenneth E. Jansen     &               ei,              h,              alfap,
15159599516SKenneth E. Jansen     &               betaT,           cp,             rk,
15259599516SKenneth E. Jansen     &               u1,              u2,             u3,
15359599516SKenneth E. Jansen     &               A0,              A1,
15459599516SKenneth E. Jansen     &               A2,              A3,
15559599516SKenneth E. Jansen     &               rLyi(:,1),       rLyi(:,2),      rLyi(:,3),  ! work arrays
15659599516SKenneth E. Jansen     &               rLyi(:,4),       rLyi(:,5),      A0DC,
15759599516SKenneth E. Jansen     &               A0inv,           dVdY)
15859599516SKenneth E. Jansen        ttim(9) = ttim(9) + tmr()
15959599516SKenneth E. Jansenc
16059599516SKenneth E. Jansenc.... calculate the convective contribution (Galerkin)
16159599516SKenneth E. Jansenc
16259599516SKenneth E. Jansen        ttim(14) = ttim(14) - secs(0.0)
16359599516SKenneth E. Jansen        call e3conv (g1yi,            g2yi,            g3yi,
16459599516SKenneth E. Jansen     &               A1,              A2,              A3,
16559599516SKenneth E. Jansen     &               rho,             pres,            T,
16659599516SKenneth E. Jansen     &               ei,              rk,              u1,
16759599516SKenneth E. Jansen     &               u2,              u3,              rLyi,
16859599516SKenneth E. Jansen     &               ri,              rmi,             EGmass,
16959599516SKenneth E. Jansen     &               shg,             shape,           WdetJ)
17059599516SKenneth E. Jansen        ttim(14) = ttim(14) + secs(0.0)
17159599516SKenneth E. Jansenc
17259599516SKenneth E. Jansenc.... calculate the diffusion contribution
17359599516SKenneth E. Jansenc
17459599516SKenneth E. Jansen        ttim(15) = ttim(15) - secs(0.0)
17559599516SKenneth E. Jansen        compK = zero
17659599516SKenneth E. Jansen        if (Navier .eq. 1) then
17759599516SKenneth E. Jansen        call e3visc (g1yi,            g2yi,            g3yi,
17859599516SKenneth E. Jansen     &               dxidx,
17959599516SKenneth E. Jansen     &               rmu,             rlm,             rlm2mu,
18059599516SKenneth E. Jansen     &               u1,              u2,              u3,
18159599516SKenneth E. Jansen     &               ri,              rmi,             stiff,
18259599516SKenneth E. Jansen     &               con, rlsli,     compK, T)
18359599516SKenneth E. Jansen        endif
18459599516SKenneth E. Jansen        ttim(15) = ttim(15) + secs(0.0)
18559599516SKenneth E. Jansenc
18659599516SKenneth E. Jansenc.... calculate the body force contribution
18759599516SKenneth E. Jansenc
18859599516SKenneth E. Jansen        if(isurf .ne. 1 .and. matflg(5,1).gt.0) then
18959599516SKenneth E. Jansen        call e3source (ri,            rmi,           rlyi,
19059599516SKenneth E. Jansen     &                 rho,           u1,            u2,
19159599516SKenneth E. Jansen     &                 u3,            pres,          sforce,
19259599516SKenneth E. Jansen     &                 dui,           dxidx,         ytargetl,
19359599516SKenneth E. Jansen     &                 xl,            shape,         bcool)
19459599516SKenneth E. Jansen        else
19559599516SKenneth E. Jansen           bcool=zero
19659599516SKenneth E. Jansen        endif
19759599516SKenneth E. Jansenc
19859599516SKenneth E. Jansenc.... calculate the least-squares contribution
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansen        ttim(16) = ttim(16) - secs(0.0)
20159599516SKenneth E. Jansen        call e3LS   (A1,              A2,            A3,
20259599516SKenneth E. Jansen     &               rho,             rmu,           cp,
20359599516SKenneth E. Jansen     &               cv,              con,           T,
20459599516SKenneth E. Jansen     &               u1,              u2,            u3,
20559599516SKenneth E. Jansen     &               rLyi,            dxidx,         tau,
20659599516SKenneth E. Jansen     &               ri,              rmi,           rk,
20759599516SKenneth E. Jansen     &               dui,             aci,           A0,
20859599516SKenneth E. Jansen     &               divqi,           shape,         shg,
20959599516SKenneth E. Jansen     &               EGmass,          stiff,         WdetJ,
21059599516SKenneth E. Jansen     &               giju,            rTLS,          raLS,
21159599516SKenneth E. Jansen     &               A0inv,           dVdY,          rerrl,
21259599516SKenneth E. Jansen     &               compK,           pres,          PTau)
21359599516SKenneth E. Jansen        ttim(16) = ttim(16) + secs(0.0)
21459599516SKenneth E. Jansenc
21559599516SKenneth E. Jansenc....  Discontinuity capturing
21659599516SKenneth E. Jansenc
21759599516SKenneth E. Jansen        if(iDC.ne.0) then
21859599516SKenneth E. Jansen          call e3dc  (g1yi,          g2yi,          g3yi,
21959599516SKenneth E. Jansen     &                A0,            raLS,          rTLS,
22059599516SKenneth E. Jansen     &                giju,          DC,
22159599516SKenneth E. Jansen     &                ri,            rmi,           stiff, A0DC)
222*f0b806cbSKenneth E. Jansen         if((intp.eq.1).and.(ierrcalc.eq.1).and.(nitr.eq.iter))  then
223*f0b806cbSKenneth E. Jansen          do i=1,npro
224*f0b806cbSKenneth E. Jansen             Tmax=maxval(yl(i,:,5))
225*f0b806cbSKenneth E. Jansen             Tmin=minval(yl(i,:,5))
226*f0b806cbSKenneth E. Jansen             rerrl(i,:,6)=(Tmax-Tmin)/T(i)
227*f0b806cbSKenneth E. Jansen          enddo
228*f0b806cbSKenneth E. Jansen!          do j=1,nshl
229*f0b806cbSKenneth E. Jansen!            rerrl(:,j,6)=rerrl(:,j,6)+DC(:,intp)  !super hack to get error indicator for shock capturing
230*f0b806cbSKenneth E. Jansen!          enddo
231*f0b806cbSKenneth E. Jansen         endif
23259599516SKenneth E. Jansen        endif
23359599516SKenneth E. Jansenc
23459599516SKenneth E. Jansenc
23559599516SKenneth E. Jansenc.... calculate the time derivative (mass) contribution to RHS
23659599516SKenneth E. Jansenc
23759599516SKenneth E. Jansen        if (ngauss .eq. 1 .and. nshl .eq. 4) then    ! trilinear tets
23859599516SKenneth E. Jansen           ttim(17) = ttim(17) - secs(0.0)
23959599516SKenneth E. Jansen           call e3juel (yl, acl, Sclr, A0, WdetJ, rl, rml)
24059599516SKenneth E. Jansen           ttim(17) = ttim(17) + secs(0.0)
24159599516SKenneth E. Jansen        else
24259599516SKenneth E. Jansen           call e3massr (aci, dui, ri,  rmi, A0)
24359599516SKenneth E. Jansen        endif
24459599516SKenneth E. Jansen
24559599516SKenneth E. Jansenc
24659599516SKenneth E. Jansenc.... calculate the time (mass) contribution to the LHS
24759599516SKenneth E. Jansenc
24859599516SKenneth E. Jansen        if (lhs .eq. 1) then
24959599516SKenneth E. Jansen           call e3massl (bcool,shape,  WdetJ,   A0,  EGmass)
25059599516SKenneth E. Jansen        endif
25159599516SKenneth E. Jansenc
25259599516SKenneth E. Jansenc....  calculate the preconditioner all at once now instead of in separate
25359599516SKenneth E. Jansenc      routines
25459599516SKenneth E. Jansenc
25559599516SKenneth E. Jansen       if(iprec.eq.1 .and. lhs.ne.1)then
25659599516SKenneth E. Jansen          ttim(18) = ttim(18) - secs(0.0)
25759599516SKenneth E. Jansen
25859599516SKenneth E. Jansen          if (itau.lt.10) then
25959599516SKenneth E. Jansen
26059599516SKenneth E. Jansen             call e3bdg(shape,       shg,             WdetJ,
26159599516SKenneth E. Jansen     &		  A1,	       A2,	        A3,
26259599516SKenneth E. Jansen     &		  A0,          bcool,           tau,
26359599516SKenneth E. Jansen     &            u1,          u2,              u3,
26459599516SKenneth E. Jansen     &            BDiagl,
26559599516SKenneth E. Jansen     &            rmu,         rlm2mu,          con)
26659599516SKenneth E. Jansen
26759599516SKenneth E. Jansen          else
26859599516SKenneth E. Jansen
26959599516SKenneth E. Jansen             call e3bdg_nd(shape,       shg,             WdetJ,
27059599516SKenneth E. Jansen     &		  A1,	       A2,	        A3,
27159599516SKenneth E. Jansen     &		  A0,          bcool,           PTau,
27259599516SKenneth E. Jansen     &            u1,          u2,              u3,
27359599516SKenneth E. Jansen     &            BDiagl,
27459599516SKenneth E. Jansen     &            rmu,         rlm2mu,          con)
27559599516SKenneth E. Jansen
27659599516SKenneth E. Jansen          endif
27759599516SKenneth E. Jansen
27859599516SKenneth E. Jansen       ttim(18) = ttim(18) + secs(0.0)
27959599516SKenneth E. Jansen       endif
28059599516SKenneth E. Jansenc
28159599516SKenneth E. Jansenc
28259599516SKenneth E. Jansenc.... multiply flux terms by shape functions and derivatives (from weight space for RHS and
28359599516SKenneth E. Jansenc     by both the weight space and solution space for the LHS stiffness term)
28459599516SKenneth E. Jansenc
28559599516SKenneth E. Jansen       ttim(19) = ttim(19) - secs(0.0)
28659599516SKenneth E. Jansen       call e3wmlt (shape,         shg,       WdetJ,
28759599516SKenneth E. Jansen     &              ri,            rmi,       rl,
28859599516SKenneth E. Jansen     &              rml,           stiff,     EGmass)
28959599516SKenneth E. Jansen       ttim(19) = ttim(19) + secs(0.0)
29059599516SKenneth E. Jansenc
29159599516SKenneth E. Jansenc.... end of integration loop
29259599516SKenneth E. Jansenc
29359599516SKenneth E. Jansen      enddo
29459599516SKenneth E. Jansen
29559599516SKenneth E. Jansen      ttim(6) = ttim(6) + secs(0.0)
29659599516SKenneth E. Jansenc
29759599516SKenneth E. Jansenc.... return
29859599516SKenneth E. Jansenc
29959599516SKenneth E. Jansen      return
30059599516SKenneth E. Jansen      end
30159599516SKenneth E. Jansenc
30259599516SKenneth E. Jansenc
30359599516SKenneth E. Jansenc
30459599516SKenneth E. Jansenc
30559599516SKenneth E. Jansen       subroutine e3Sclr (ycl,          acl,
30659599516SKenneth E. Jansen     &                    dwl,         elDwl,
30759599516SKenneth E. Jansen     &                    shp,         sgn,
30859599516SKenneth E. Jansen     &                    shgl,        xl,
30959599516SKenneth E. Jansen     &                    rtl,         rmtl,
31059599516SKenneth E. Jansen     &                    qtl,         EGmasst)
31159599516SKenneth E. Jansenc
31259599516SKenneth E. Jansenc----------------------------------------------------------------------
31359599516SKenneth E. Jansenc
31459599516SKenneth E. Jansenc This routine is the 3D element routine for the N-S equations.
31559599516SKenneth E. Jansenc This routine calculates the RHS residual and if requested the
31659599516SKenneth E. Jansenc modified residual or the LHS consistent mass matrix or the block-
31759599516SKenneth E. Jansenc diagonal preconditioner.
31859599516SKenneth E. Jansenc
31959599516SKenneth E. Jansenc input:    e    a   1..5   when we think of U^e_a  and U is 5 variables
32059599516SKenneth E. Jansenc  ycl      (npro,nshl,ndof)     : Y variables
32159599516SKenneth E. Jansenc  actl    (npro,nshl)          : scalar variable time derivative
32259599516SKenneth E. Jansenc  dwl     (npro,nenl)          : distances to wall
32359599516SKenneth E. Jansenc  shp     (nen,ngauss)          : element shape-functions  N_a
32459599516SKenneth E. Jansenc  shgl    (nsd,nen,ngauss)      : element local-shape-functions N_{a,xi}
32559599516SKenneth E. Jansenc  xl      (npro,nenl,nsd)      : nodal coordinates at current step (x^e_a)
32659599516SKenneth E. Jansenc  qtl     (npro,nshl)          : diffusive flux vector (don't worry)
32759599516SKenneth E. Jansenc
32859599516SKenneth E. Jansenc output:
32959599516SKenneth E. Jansenc  rtl     (npro,nshl)          : element RHS residual    (G^e_a)
33059599516SKenneth E. Jansenc  rmtl    (npro,nshl)          : element modified residual  (G^e_a tilde)
33159599516SKenneth E. Jansenc  EGmasst (npro,nshape,nshape) : element LHS tangent mass matrix (dG^e_a
33259599516SKenneth E. Jansenc                                                                  dY_b  )
33359599516SKenneth E. Jansenc
33459599516SKenneth E. Jansenc
33559599516SKenneth E. Jansenc Note: This routine will calculate the element matrices for the
33659599516SKenneth E. Jansenc        Hulbert's generalized alpha method integrator
33759599516SKenneth E. Jansenc
33859599516SKenneth E. Jansenc Note: nedof = nflow*nshl is the total number of degrees of freedom
33959599516SKenneth E. Jansenc       at each element node.
34059599516SKenneth E. Jansenc
34159599516SKenneth E. Jansenc Mathematics done by:  Michel Mallet, Farzin Shakib (version 1)
34259599516SKenneth E. Jansenc                       Farzin Shakib                (version 2)
34359599516SKenneth E. Jansenc
34459599516SKenneth E. Jansenc
34559599516SKenneth E. Jansenc Zdenek Johan, Summer 1990.   (Modified from e2.f)
34659599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.   (Fortran 90)
34759599516SKenneth E. Jansenc Kenneth Jansen, Winter 1997. (Primitive Variables)
34859599516SKenneth E. Jansenc Chris Whiting, Winter 1998.  (LHS matrix formation)
34959599516SKenneth E. Jansenc----------------------------------------------------------------------
35059599516SKenneth E. Jansenc
35159599516SKenneth E. Jansen        include "common.h"
35259599516SKenneth E. Jansenc
35359599516SKenneth E. Jansen        dimension ycl(npro,nshl,ndof),
35459599516SKenneth E. Jansen     &            acl(npro,nshl,ndof),
35559599516SKenneth E. Jansen     &            dwl(npro,nenl),
35659599516SKenneth E. Jansen     &            shp(nshl,ngauss),         shgl(nsd,nshl,ngauss),
35759599516SKenneth E. Jansen     &            xl(npro,nenl,nsd),
35859599516SKenneth E. Jansen     &            rtl(npro,nshl),           qtl(npro,nshl),
35959599516SKenneth E. Jansen     &            rmtl(npro,nshl),          Diagl(npro,nshl),
36059599516SKenneth E. Jansen     &            EGmasst(npro,nshape,nshape),
36159599516SKenneth E. Jansen     &            dist2w(npro),             sgn(npro,nshl),
362513954efSKenneth E. Jansen     &            vort(npro),               gVnrm(npro),
36359599516SKenneth E. Jansen     &            rmu(npro),                con(npro),
36459599516SKenneth E. Jansen     &            T(npro),                  cp(npro),
36559599516SKenneth E. Jansen     &            g1yti(npro),              acti(npro),
36659599516SKenneth E. Jansen     &            g2yti(npro),              g3yti(npro),
36759599516SKenneth E. Jansen     &            Sclr(npro),               srcp (npro)
36859599516SKenneth E. Jansen
36959599516SKenneth E. Jansenc
37059599516SKenneth E. Jansen        dimension shg(npro,nshl,nsd)
37159599516SKenneth E. Jansenc
37259599516SKenneth E. Jansen        dimension dxidx(npro,nsd,nsd),     WdetJ(npro)
37359599516SKenneth E. Jansenc
37459599516SKenneth E. Jansen        dimension rho(npro),               rk(npro),
37559599516SKenneth E. Jansen     &            u1(npro),                u2(npro),
37659599516SKenneth E. Jansen     &            u3(npro)
37759599516SKenneth E. Jansenc
37859599516SKenneth E. Jansen        dimension A0t(npro),               A1t(npro),
37959599516SKenneth E. Jansen     &            A2t(npro),               A3t(npro)
38059599516SKenneth E. Jansenc
38159599516SKenneth E. Jansen        dimension rLyti(npro),             raLSt(npro),
38259599516SKenneth E. Jansen     &            rTLSt(npro),             giju(npro,6),
38359599516SKenneth E. Jansen     &            DCt(npro, ngauss)
38459599516SKenneth E. Jansenc
38559599516SKenneth E. Jansen        dimension rti(npro,nsd+1),         rmti(npro,nsd+1),
38659599516SKenneth E. Jansen     &            stifft(npro,nsd,nsd),
38759599516SKenneth E. Jansen     &            shape(npro,nshl),        shdrv(npro,nsd,nshl)
38859599516SKenneth E. Jansen        real*8    elDwl(npro)
38959599516SKenneth E. Jansenc
39059599516SKenneth E. Jansen        ttim(6) = ttim(6) - tmr()
39159599516SKenneth E. Jansenc
39259599516SKenneth E. Jansenc.... loop through the integration points
39359599516SKenneth E. Jansenc
39459599516SKenneth E. Jansen        elDwl(:)=zero
39559599516SKenneth E. Jansen        do intp = 1, ngauss
39659599516SKenneth E. Jansenc
39759599516SKenneth E. Jansenc.... if Det. .eq. 0, do not include this point
39859599516SKenneth E. Jansenc
39959599516SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
40059599516SKenneth E. Jansenc
40159599516SKenneth E. Jansenc
40259599516SKenneth E. Jansenc.... create a matrix of shape functions (and derivatives) for each
40359599516SKenneth E. Jansenc     element at this quadrature point. These arrays will contain
40459599516SKenneth E. Jansenc     the correct signs for the hierarchic basis
40559599516SKenneth E. Jansenc
40659599516SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
40759599516SKenneth E. Jansen     &              shape,        shdrv)
40859599516SKenneth E. Jansenc
40959599516SKenneth E. Jansenc.... initialize
41059599516SKenneth E. Jansenc
41159599516SKenneth E. Jansen        rlyti = zero
41259599516SKenneth E. Jansen        rti   = zero
41359599516SKenneth E. Jansen        rmti  = zero
41459599516SKenneth E. Jansen        srcp  = zero
41559599516SKenneth E. Jansen        stifft = zero
41659599516SKenneth E. Jansenc        if (lhs .eq. 1) stifft = zero
41759599516SKenneth E. Jansenc
41859599516SKenneth E. Jansenc
41959599516SKenneth E. Jansenc.... calculate the integration variables
42059599516SKenneth E. Jansenc
42159599516SKenneth E. Jansen        ttim(8) = ttim(8) - tmr()
42259599516SKenneth E. Jansenc
42359599516SKenneth E. Jansen        call e3ivarSclr(ycl,              acl,        acti,
42459599516SKenneth E. Jansen     &                  shape,           shdrv,      xl,
42559599516SKenneth E. Jansen     &                  T,               cp,
42659599516SKenneth E. Jansen     &                  dxidx,           Sclr,
427513954efSKenneth E. Jansen     &                  Wdetj,           vort,       gVnrm,
42859599516SKenneth E. Jansen     &                  g1yti,           g2yti,      g3yti,
42959599516SKenneth E. Jansen     &                  rho,             rmu,        con,
43059599516SKenneth E. Jansen     &                  rk,              u1,         u2,
43159599516SKenneth E. Jansen     &                  u3,              shg,        dwl,
43259599516SKenneth E. Jansen     &                  dist2w)
43359599516SKenneth E. Jansen        ttim(8) = ttim(8) + tmr()
43459599516SKenneth E. Jansen
43559599516SKenneth E. Jansenc
43659599516SKenneth E. Jansenc.... calculate the source term contribution
43759599516SKenneth E. Jansenc
43859599516SKenneth E. Jansen        if(nosource.ne.1)
43959599516SKenneth E. Jansen     &  call e3sourceSclr (Sclr,    rho,    rmu,
440513954efSKenneth E. Jansen     &                     dist2w,  vort,   gVnrm,   con,
44159599516SKenneth E. Jansen     &                     g1yti,  g2yti,   g3yti,
44259599516SKenneth E. Jansen     &                     rti,    rLyti,   srcp,
44359599516SKenneth E. Jansen     &                     ycl,     shape,   u1,
44459599516SKenneth E. Jansen     &                     u2,     u3,	     xl,
44559599516SKenneth E. Jansen     &                     elDwl)
44659599516SKenneth E. Jansenc
44759599516SKenneth E. Jansen         if (ilset.eq.2 .and. isclr.eq.2) then
44859599516SKenneth E. Jansen          rk = pt5 * ( u1**2 + u2**2 + u3**2 )
44959599516SKenneth E. Jansen         endif
45059599516SKenneth E. Jansenc.... calculate the relevant matrices
45159599516SKenneth E. Jansenc
45259599516SKenneth E. Jansen        ttim(9) = ttim(9) - tmr()
45359599516SKenneth E. Jansen        call e3mtrxSclr (rho,
45459599516SKenneth E. Jansen     &                   u1,            u2,         u3,
45559599516SKenneth E. Jansen     &                   A0t,           A1t,
45659599516SKenneth E. Jansen     &                   A2t,           A3t)
45759599516SKenneth E. Jansen        ttim(9) = ttim(9) + tmr()
45859599516SKenneth E. Jansenc
45959599516SKenneth E. Jansenc.... calculate the convective contribution (Galerkin)
46059599516SKenneth E. Jansenc
46159599516SKenneth E. Jansen        ttim(14) = ttim(14) - tmr()
46259599516SKenneth E. Jansen        call e3convSclr (g1yti,        g2yti,       g3yti,
46359599516SKenneth E. Jansen     &                   A1t,          A2t,         A3t,
46459599516SKenneth E. Jansen     &                   rho,          u1,          Sclr,
46559599516SKenneth E. Jansen     &                   u2,           u3,          rLyti,
46659599516SKenneth E. Jansen     &                   rti,          rmti,        EGmasst,
46759599516SKenneth E. Jansen     &                   shg,          shape,       WdetJ)
46859599516SKenneth E. Jansen        ttim(14) = ttim(14) + tmr()
46959599516SKenneth E. Jansenc
47059599516SKenneth E. Jansenc.... calculate the diffusion contribution
47159599516SKenneth E. Jansenc
47259599516SKenneth E. Jansen        ttim(15) = ttim(15) - tmr()
47359599516SKenneth E. Jansen        if (Navier .eq. 1) then
47459599516SKenneth E. Jansen        call e3viscSclr (g1yti,        g2yti,         g3yti,
47559599516SKenneth E. Jansen     &                   rmu,          Sclr,          rho,
47659599516SKenneth E. Jansen     &                   rti,          rmti,          stifft )
47759599516SKenneth E. Jansen        endif
47859599516SKenneth E. Jansen         ttim(15) = ttim(15) + tmr()
47959599516SKenneth E. Jansenc
48059599516SKenneth E. Jansen         if (ilset.eq.2)  srcp = zero
48159599516SKenneth E. Jansen
48259599516SKenneth E. Jansenc
48359599516SKenneth E. Jansenc.... calculate the least-squares contribution
48459599516SKenneth E. Jansenc
48559599516SKenneth E. Jansen        ttim(16) = ttim(16) - tmr()
48659599516SKenneth E. Jansen        call e3LSSclr(A1t,             A2t,             A3t,
48759599516SKenneth E. Jansen     &                rho,             rmu,             rtLSt,
48859599516SKenneth E. Jansen     &                u1,              u2,              u3,
48959599516SKenneth E. Jansen     &                rLyti,           dxidx,           raLSt,
49059599516SKenneth E. Jansen     &                rti,             rk,              giju,
49159599516SKenneth E. Jansen     &                acti,            A0t,
49259599516SKenneth E. Jansen     &                shape,           shg,
49359599516SKenneth E. Jansen     &                EGmasst,         stifft,          WdetJ,
49459599516SKenneth E. Jansen     &                srcp)
49559599516SKenneth E. Jansen        ttim(16) = ttim(16) + tmr()
49659599516SKenneth E. Jansenc
49759599516SKenneth E. Jansenc******************************DC TERMS*****************************
49859599516SKenneth E. Jansen        if (idcsclr(1) .ne. 0) then
49959599516SKenneth E. Jansen           if ((idcsclr(2).eq.1 .and. isclr.eq.1) .or.
50059599516SKenneth E. Jansen     &         (idcsclr(2).eq.2 .and. isclr.eq.2)) then   ! scalar with dc
50159599516SKenneth E. Jansen              call e3dcSclr(g1yti, g2yti,          g3yti,
50259599516SKenneth E. Jansen     &             A0t,            raLSt,          rTLSt,
50359599516SKenneth E. Jansen     &             DCt,            giju,
50459599516SKenneth E. Jansen     &             rti,            rmti,           stifft)
50559599516SKenneth E. Jansen           endif
50659599516SKenneth E. Jansen        endif
50759599516SKenneth E. Jansenc
50859599516SKenneth E. Jansenc******************************************************************
50959599516SKenneth E. Jansenc.... calculate the time derivative (mass) contribution to RHS
51059599516SKenneth E. Jansenc
51159599516SKenneth E. Jansen
51259599516SKenneth E. Jansen           call e3massrSclr (acti, rti, A0t)
51359599516SKenneth E. Jansenc
51459599516SKenneth E. Jansenc.... calculate the time (mass) contribution to the LHS
51559599516SKenneth E. Jansenc
51659599516SKenneth E. Jansen         if (lhs .eq. 1) then
51759599516SKenneth E. Jansen            call e3masslSclr (shape,  WdetJ,   A0t,  EGmasst,srcp)
51859599516SKenneth E. Jansen         endif
51959599516SKenneth E. Jansenc
52059599516SKenneth E. Jansen
52159599516SKenneth E. Jansenc.... multiply flux terms by shape functions and derivatives (from weight space for RHS and
52259599516SKenneth E. Jansenc     by both the weight space and solution space for the LHS stiffness term)
52359599516SKenneth E. Jansenc
52459599516SKenneth E. Jansen       ttim(19) = ttim(19) - tmr()
52559599516SKenneth E. Jansen       call e3wmltSclr (shape,           shg,             WdetJ,
52659599516SKenneth E. Jansen     &                  rti,
52759599516SKenneth E. Jansen     &                  rtl,             stifft,          EGmasst)
52859599516SKenneth E. Jansen       ttim(19) = ttim(19) + tmr()
52959599516SKenneth E. Jansenc
53059599516SKenneth E. Jansenc.... end of the loop
53159599516SKenneth E. Jansenc
53259599516SKenneth E. Jansen      enddo
53359599516SKenneth E. Jansen	qptinv=one/ngauss
53459599516SKenneth E. Jansen	elDwl(:)=elDwl(:)*qptinv
53559599516SKenneth E. Jansen
53659599516SKenneth E. Jansen
53759599516SKenneth E. Jansen      ttim(6) = ttim(6) + tmr()
53859599516SKenneth E. Jansenc
53959599516SKenneth E. Jansenc.... return
54059599516SKenneth E. Jansenc
54159599516SKenneth E. Jansen      return
54259599516SKenneth E. Jansen      end
54359599516SKenneth E. Jansen
544