xref: /phasta/phSolver/compressible/e3ivar.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
159599516SKenneth E. Jansen        subroutine e3ivar (yl,      ycl,     acl,
259599516SKenneth E. Jansen     &                     Sclr,    shape,   shgl,
359599516SKenneth E. Jansen     &                     xl,      dui,     aci,
459599516SKenneth E. Jansen     &                     g1yi,    g2yi,    g3yi,
559599516SKenneth E. Jansen     &                     shg,     dxidx,   WdetJ,
659599516SKenneth E. Jansen     &                     rho,     pres,    T,
759599516SKenneth E. Jansen     &                     ei,      h,       alfap,
859599516SKenneth E. Jansen     &                     betaT,   cp,      rk,
959599516SKenneth E. Jansen     &                     u1,      u2,      u3,
1059599516SKenneth E. Jansen     &                     ql,      divqi,   sgn, tmp,
1159599516SKenneth E. Jansen     &                     rmu,     rlm,     rlm2mu,
1259599516SKenneth E. Jansen     &                     con,     rlsl,    rlsli,
1359599516SKenneth E. Jansen     &                     xmudmi,  sforce,  cv)
1459599516SKenneth E. Jansenc
1559599516SKenneth E. Jansenc----------------------------------------------------------------------
1659599516SKenneth E. Jansenc
1759599516SKenneth E. Jansenc  This routine computes the variables at integration point.
1859599516SKenneth E. Jansenc
1959599516SKenneth E. Jansenc input:
2059599516SKenneth E. Jansenc  yl     (npro,nshl,nflow)     : primitive variables (NO SCALARS)
2159599516SKenneth E. Jansenc  ycl    (npro,nshl,ndof)      : primitive variables at current step
2259599516SKenneth E. Jansenc  acl    (npro,nshl,ndof)      : prim.var. accel. at the current step
2359599516SKenneth E. Jansenc  shape  (npro,nshl)           : element shape-functions
2459599516SKenneth E. Jansenc  shgl   (nsd,nen)             : element local-grad-shape-functions
2559599516SKenneth E. Jansenc  xl     (npro,nenl,nsd)       : nodal coordinates at current step
2659599516SKenneth E. Jansenc  ql     (npro,nshl,(nflow-1)*nsd) : diffusive flux vector
2759599516SKenneth E. Jansenc  rlsl   (npro,nshl,6)       : resolved Leonard stresses
2859599516SKenneth E. Jansenc  sgn    (npro,nshl)         : signs of shape functions
2959599516SKenneth E. Jansenc
3059599516SKenneth E. Jansenc output:
3159599516SKenneth E. Jansenc  dui    (npro,nflow)           : delta U variables at current step
3259599516SKenneth E. Jansenc  aci    (npro,nflow)           : primvar accel. variables at current step
3359599516SKenneth E. Jansenc  g1yi   (npro,nflow)           : grad-y in direction 1
3459599516SKenneth E. Jansenc  g2yi   (npro,nflow)           : grad-y in direction 2
3559599516SKenneth E. Jansenc  g3yi   (npro,nflow)           : grad-y in direction 3
3659599516SKenneth E. Jansenc  shg    (npro,nshl,nsd)       : element global grad-shape-functions
3759599516SKenneth E. Jansenc  dxidx  (npro,nsd,nsd)        : inverse of deformation gradient
3859599516SKenneth E. Jansenc  WdetJ  (npro)                : weighted Jacobian
3959599516SKenneth E. Jansenc  rho    (npro)                : density
4059599516SKenneth E. Jansenc  pres   (npro)                : pressure
4159599516SKenneth E. Jansenc  T      (npro)                : temperature
4259599516SKenneth E. Jansenc  ei     (npro)                : internal energy
4359599516SKenneth E. Jansenc  h      (npro)                : enthalpy
4459599516SKenneth E. Jansenc  alfap  (npro)                : expansivity
4559599516SKenneth E. Jansenc  betaT  (npro)                : isothermal compressibility
4659599516SKenneth E. Jansenc  cp     (npro)                : specific heat at constant pressure
4759599516SKenneth E. Jansenc  rk     (npro)                : kinetic energy
4859599516SKenneth E. Jansenc  u1     (npro)                : x1-velocity component
4959599516SKenneth E. Jansenc  u2     (npro)                : x2-velocity component
5059599516SKenneth E. Jansenc  u3     (npro)                : x3-velocity component
5159599516SKenneth E. Jansenc  divqi  (npro,nflow-1)        : divergence of diffusive flux
5259599516SKenneth E. Jansenc  rlsli  (npro,6)              : resolved Leonard stresses at quad pt
5359599516SKenneth E. Jansenc
5459599516SKenneth E. Jansenc Zdenek Johan, Summer 1990. (Modified from e2ivar.f)
5559599516SKenneth E. Jansenc Zdenek Johan, Winter 1991. (Fortran 90)
5659599516SKenneth E. Jansenc Kenneth Jansen, Winter 1997. Primitive Variables
5759599516SKenneth E. Jansenc----------------------------------------------------------------------
5859599516SKenneth E. Jansenc
5959599516SKenneth E. Jansen        include "common.h"
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansenc  passed arrays
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansen        dimension yl(npro,nshl,nflow),        ycl(npro,nshl,ndof),
6459599516SKenneth E. Jansen     &            acl(npro,nshl,ndof),
6559599516SKenneth E. Jansen     &            shape(npro,nshl),
6659599516SKenneth E. Jansen     &            shgl(npro,nsd,nshl),      xl(npro,nenl,nsd),
6759599516SKenneth E. Jansen     &            dui(npro,nflow),
6859599516SKenneth E. Jansen     &            aci(npro,nflow),            g1yi(npro,nflow),
6959599516SKenneth E. Jansen     &            g2yi(npro,nflow),           g3yi(npro,nflow),
7059599516SKenneth E. Jansen     &            shg(npro,nshl,nsd),        dxidx(npro,nsd,nsd),
7159599516SKenneth E. Jansen     &            WdetJ(npro),
7259599516SKenneth E. Jansen     &            rho(npro),                 pres(npro),
7359599516SKenneth E. Jansen     &            T(npro),                   ei(npro),
7459599516SKenneth E. Jansen     &            h(npro),                   alfap(npro),
7559599516SKenneth E. Jansen     &            betaT(npro),               cp(npro),
7659599516SKenneth E. Jansen     &            rk(npro),                  cv(npro),
7759599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
7859599516SKenneth E. Jansen     &            u3(npro),                  divqi(npro,nflow),
7959599516SKenneth E. Jansen     &            ql(npro,nshl,idflx),
8059599516SKenneth E. Jansen     &            rmu(npro),                 rlm(npro),
8159599516SKenneth E. Jansen     &            rlm2mu(npro),              con(npro),
8259599516SKenneth E. Jansen     &            Sclr(npro)
8359599516SKenneth E. Jansenc
8459599516SKenneth E. Jansen        dimension tmp(npro),  dxdxi(npro,nsd,nsd),  sgn(npro,nshl)
8559599516SKenneth E. Jansen        dimension rlsl(npro,nshl,6),         rlsli(npro,6),
8659599516SKenneth E. Jansen     &            xmudmi(npro,ngauss)
8759599516SKenneth E. Jansen        dimension gyti(npro,nsd),            gradh(npro,nsd),
8859599516SKenneth E. Jansen     &            sforce(npro,3),            weber(npro)
8959599516SKenneth E. Jansen        ttim(20) = ttim(20) - secs(0.0)
9059599516SKenneth E. Jansen
9159599516SKenneth E. Jansenc
9259599516SKenneth E. Jansen        ttim(10) = ttim(10) - secs(0.0)
9359599516SKenneth E. Jansen
9459599516SKenneth E. Jansen        dui = zero
9559599516SKenneth E. Jansenc
9659599516SKenneth E. Jansen        do n = 1, nshl
9759599516SKenneth E. Jansen           dui(:,1) = dui(:,1) + shape(:,n) * yl(:,n,1) ! p
9859599516SKenneth E. Jansen           dui(:,2) = dui(:,2) + shape(:,n) * yl(:,n,2) ! u1
9959599516SKenneth E. Jansen           dui(:,3) = dui(:,3) + shape(:,n) * yl(:,n,3) ! u2
10059599516SKenneth E. Jansen           dui(:,4) = dui(:,4) + shape(:,n) * yl(:,n,4) ! u3
10159599516SKenneth E. Jansen           dui(:,5) = dui(:,5) + shape(:,n) * yl(:,n,5) ! T
10259599516SKenneth E. Jansen        enddo
10359599516SKenneth E. Jansenc
104*f4d0b58bSKenneth E. Jansen!      flops = flops + 10*nshl*npro
10559599516SKenneth E. Jansenc
10659599516SKenneth E. Jansenc
10759599516SKenneth E. Jansenc.... compute conservative variables
10859599516SKenneth E. Jansenc
10959599516SKenneth E. Jansen        rk = pt5 * (dui(:,2)**2 + dui(:,3)**2 + dui(:,4)**2)
11059599516SKenneth E. Jansenc
11159599516SKenneth E. Jansen        if (iLSet .ne. 0)then
11259599516SKenneth E. Jansen           Sclr = zero
11359599516SKenneth E. Jansen           isc=abs(iRANS)+6
11459599516SKenneth E. Jansen           do n = 1, nshl
11559599516SKenneth E. Jansen              Sclr = Sclr + shape(:,n) * ycl(:,n,isc)
11659599516SKenneth E. Jansen           enddo
11759599516SKenneth E. Jansen        endif
11859599516SKenneth E. Jansen
11959599516SKenneth E. Jansen        ithm = 6
12059599516SKenneth E. Jansen        call getthm (dui(:,1),   dui(:,5),     Sclr,
12159599516SKenneth E. Jansen     &               rk,         rho,          ei,
12259599516SKenneth E. Jansen     &               tmp,        tmp,          tmp,
12359599516SKenneth E. Jansen     &               tmp,        tmp,          tmp,
12459599516SKenneth E. Jansen     &               tmp,        tmp)
12559599516SKenneth E. Jansenc
12659599516SKenneth E. Jansenc
12759599516SKenneth E. Jansen        dui(:,1) = rho
12859599516SKenneth E. Jansen        dui(:,2) = rho * dui(:,2)
12959599516SKenneth E. Jansen        dui(:,3) = rho * dui(:,3)
13059599516SKenneth E. Jansen        dui(:,4) = rho * dui(:,4)
13159599516SKenneth E. Jansen        dui(:,5) = rho * (ei + rk)
13259599516SKenneth E. Jansen
13359599516SKenneth E. Jansen
13459599516SKenneth E. Jansen       ttim(10) = ttim(10) + secs(0.0)
13559599516SKenneth E. Jansenc
13659599516SKenneth E. Jansenc.... ------------->  Primitive variables at int. point  <--------------
13759599516SKenneth E. Jansenc
13859599516SKenneth E. Jansenc.... compute primitive variables
13959599516SKenneth E. Jansenc
14059599516SKenneth E. Jansen       ttim(11) = ttim(11) - secs(0.0)
14159599516SKenneth E. Jansen
14259599516SKenneth E. Jansen       pres = zero
14359599516SKenneth E. Jansen       u1   = zero
14459599516SKenneth E. Jansen       u2   = zero
14559599516SKenneth E. Jansen       u3   = zero
14659599516SKenneth E. Jansen       T    = zero
14759599516SKenneth E. Jansen       do n = 1, nshl
14859599516SKenneth E. Jansenc
14959599516SKenneth E. Jansenc  y(int)=SUM_{a=1}^nshl (N_a(int) Ya)
15059599516SKenneth E. Jansenc
15159599516SKenneth E. Jansen          pres = pres + shape(:,n) * ycl(:,n,1)
15259599516SKenneth E. Jansen          u1   = u1   + shape(:,n) * ycl(:,n,2)
15359599516SKenneth E. Jansen          u2   = u2   + shape(:,n) * ycl(:,n,3)
15459599516SKenneth E. Jansen          u3   = u3   + shape(:,n) * ycl(:,n,4)
15559599516SKenneth E. Jansen          T    = T    + shape(:,n) * ycl(:,n,5)
15659599516SKenneth E. Jansen       enddo
15759599516SKenneth E. Jansen
15859599516SKenneth E. Jansen       if( (iLES.gt.10).and.(iLES.lt.20))  then  ! bardina
15959599516SKenneth E. Jansen       rlsli = zero
16059599516SKenneth E. Jansen       do n = 1, nshl
16159599516SKenneth E. Jansen
16259599516SKenneth E. Jansen          rlsli(:,1) = rlsli(:,1) + shape(:,n) * rlsl(:,n,1)
16359599516SKenneth E. Jansen          rlsli(:,2) = rlsli(:,2) + shape(:,n) * rlsl(:,n,2)
16459599516SKenneth E. Jansen          rlsli(:,3) = rlsli(:,3) + shape(:,n) * rlsl(:,n,3)
16559599516SKenneth E. Jansen          rlsli(:,4) = rlsli(:,4) + shape(:,n) * rlsl(:,n,4)
16659599516SKenneth E. Jansen          rlsli(:,5) = rlsli(:,5) + shape(:,n) * rlsl(:,n,5)
16759599516SKenneth E. Jansen          rlsli(:,6) = rlsli(:,6) + shape(:,n) * rlsl(:,n,6)
16859599516SKenneth E. Jansen
16959599516SKenneth E. Jansen       enddo
17059599516SKenneth E. Jansen       else
17159599516SKenneth E. Jansen          rlsli = zero
17259599516SKenneth E. Jansen       endif
17359599516SKenneth E. Jansenc
17459599516SKenneth E. Jansen
17559599516SKenneth E. Jansen       ttim(11) = ttim(11) + secs(0.0)
17659599516SKenneth E. Jansen
17759599516SKenneth E. Jansenc
17859599516SKenneth E. Jansenc.... ----------------------->  accel. at int. point  <----------------------
17959599516SKenneth E. Jansenc
18059599516SKenneth E. Jansen
18159599516SKenneth E. Jansenc       if (ires .ne. 2)  then
18259599516SKenneth E. Jansen          ttim(12) = ttim(12) - secs(0.0)
18359599516SKenneth E. Jansenc
18459599516SKenneth E. Jansenc.... compute primitive variables
18559599516SKenneth E. Jansenc
18659599516SKenneth E. Jansen          aci = zero
18759599516SKenneth E. Jansenc
18859599516SKenneth E. Jansen          do n = 1, nshl
18959599516SKenneth E. Jansen             aci(:,1) = aci(:,1) + shape(:,n) * acl(:,n,1)
19059599516SKenneth E. Jansen             aci(:,2) = aci(:,2) + shape(:,n) * acl(:,n,2)
19159599516SKenneth E. Jansen             aci(:,3) = aci(:,3) + shape(:,n) * acl(:,n,3)
19259599516SKenneth E. Jansen             aci(:,4) = aci(:,4) + shape(:,n) * acl(:,n,4)
19359599516SKenneth E. Jansen             aci(:,5) = aci(:,5) + shape(:,n) * acl(:,n,5)
19459599516SKenneth E. Jansen          enddo
19559599516SKenneth E. Jansenc
196*f4d0b58bSKenneth E. Jansen!      flops = flops + 10*nshl*npro
19759599516SKenneth E. Jansen          ttim(12) = ttim(12) + secs(0.0)
19859599516SKenneth E. Jansenc       endif                    !ires .ne. 2
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansenc.... ----------------->  Thermodynamic Properties  <-------------------
20159599516SKenneth E. Jansenc
20259599516SKenneth E. Jansenc.... compute the kinetic energy
20359599516SKenneth E. Jansenc
20459599516SKenneth E. Jansen        ttim(27) = ttim(27) - secs(0.0)
20559599516SKenneth E. Jansenc
20659599516SKenneth E. Jansen        rk = pt5 * ( u1**2 + u2**2 + u3**2 )
20759599516SKenneth E. Jansenc
20859599516SKenneth E. Jansenc.... get the thermodynamic properties
20959599516SKenneth E. Jansenc
21059599516SKenneth E. Jansen        if (iLSet .ne. 0)then
21159599516SKenneth E. Jansen           Sclr = zero
21259599516SKenneth E. Jansen           isc=abs(iRANS)+6
21359599516SKenneth E. Jansen           do n = 1, nshl
21459599516SKenneth E. Jansen              Sclr = Sclr + shape(:,n) * ycl(:,n,isc)
21559599516SKenneth E. Jansen           enddo
21659599516SKenneth E. Jansen        endif
21759599516SKenneth E. Jansen
21859599516SKenneth E. Jansen        ithm = 7
21959599516SKenneth E. Jansen        call getthm (pres,            T,               Sclr,
22059599516SKenneth E. Jansen     &               rk,              rho,             ei,
22159599516SKenneth E. Jansen     &               h,               tmp,             cv,
22259599516SKenneth E. Jansen     &               cp,              alfap,           betaT,
22359599516SKenneth E. Jansen     &               tmp,             tmp)
22459599516SKenneth E. Jansenc
22559599516SKenneth E. Jansen        ttim(27) = ttim(27) + secs(0.0)
22659599516SKenneth E. Jansenc
22759599516SKenneth E. Jansenc ........Get the material properties
22859599516SKenneth E. Jansenc
22959599516SKenneth E. Jansen        call getDiff (T,        cp,       rho,        ycl,
23059599516SKenneth E. Jansen     &                rmu,      rlm,      rlm2mu,     con,  shape,
23159599516SKenneth E. Jansen     &                xmudmi,   xl)
23259599516SKenneth E. Jansen
23359599516SKenneth E. Jansenc.... --------------------->  Element Metrics  <-----------------------
23459599516SKenneth E. Jansenc
23559599516SKenneth E. Jansen      ttim(26) = ttim(26) - secs(0.0)
23659599516SKenneth E. Jansenc
23759599516SKenneth E. Jansen        call e3metric( xl,         shgl,        dxidx,
23859599516SKenneth E. Jansen     &                 shg,        WdetJ)
23959599516SKenneth E. Jansenc
24059599516SKenneth E. Jansenc
24159599516SKenneth E. Jansen        ttim(26) = ttim(26) + secs(0.0)
24259599516SKenneth E. Jansenc
24359599516SKenneth E. Jansenc.... --------------------->  Global Gradients  <-----------------------
24459599516SKenneth E. Jansenc
24559599516SKenneth E. Jansen        ttim(13) = ttim(13) - secs(0.0)
24659599516SKenneth E. Jansen
24759599516SKenneth E. Jansen        g1yi = zero
24859599516SKenneth E. Jansen        g2yi = zero
24959599516SKenneth E. Jansen        g3yi = zero
25059599516SKenneth E. Jansenc
25159599516SKenneth E. Jansen        ttim(13) = ttim(13) + secs(0.0)
25259599516SKenneth E. Jansen        ttim(7) = ttim(7) - secs(0.0)
25359599516SKenneth E. Jansenc
25459599516SKenneth E. Jansenc.... compute the global gradient of Y-variables
25559599516SKenneth E. Jansenc
25659599516SKenneth E. Jansenc
25759599516SKenneth E. Jansenc  Y_{,x_i}=SUM_{a=1}^nshl (N_{a,x_i}(int) Ya)
25859599516SKenneth E. Jansenc
25959599516SKenneth E. Jansen        if(nshl.eq.4) then
26059599516SKenneth E. Jansen          g1yi(:,1) = g1yi(:,1) + shg(:,1,1) * yl(:,1,1)
26159599516SKenneth E. Jansen     &                          + shg(:,2,1) * yl(:,2,1)
26259599516SKenneth E. Jansen     &                          + shg(:,3,1) * yl(:,3,1)
26359599516SKenneth E. Jansen     &                          + shg(:,4,1) * yl(:,4,1)
26459599516SKenneth E. Jansenc
26559599516SKenneth E. Jansen          g1yi(:,2) = g1yi(:,2) + shg(:,1,1) * yl(:,1,2)
26659599516SKenneth E. Jansen     &                          + shg(:,2,1) * yl(:,2,2)
26759599516SKenneth E. Jansen     &                          + shg(:,3,1) * yl(:,3,2)
26859599516SKenneth E. Jansen     &                          + shg(:,4,1) * yl(:,4,2)
26959599516SKenneth E. Jansenc
27059599516SKenneth E. Jansen          g1yi(:,3) = g1yi(:,3) + shg(:,1,1) * yl(:,1,3)
27159599516SKenneth E. Jansen     &                          + shg(:,2,1) * yl(:,2,3)
27259599516SKenneth E. Jansen     &                          + shg(:,3,1) * yl(:,3,3)
27359599516SKenneth E. Jansen     &                          + shg(:,4,1) * yl(:,4,3)
27459599516SKenneth E. Jansenc
27559599516SKenneth E. Jansen          g1yi(:,4) = g1yi(:,4) + shg(:,1,1) * yl(:,1,4)
27659599516SKenneth E. Jansen     &                          + shg(:,2,1) * yl(:,2,4)
27759599516SKenneth E. Jansen     &                          + shg(:,3,1) * yl(:,3,4)
27859599516SKenneth E. Jansen     &                          + shg(:,4,1) * yl(:,4,4)
27959599516SKenneth E. Jansenc
28059599516SKenneth E. Jansen          g1yi(:,5) = g1yi(:,5) + shg(:,1,1) * yl(:,1,5)
28159599516SKenneth E. Jansen     &                          + shg(:,2,1) * yl(:,2,5)
28259599516SKenneth E. Jansen     &                          + shg(:,3,1) * yl(:,3,5)
28359599516SKenneth E. Jansen     &                          + shg(:,4,1) * yl(:,4,5)
28459599516SKenneth E. Jansenc
28559599516SKenneth E. Jansen          g2yi(:,1) = g2yi(:,1) + shg(:,1,2) * yl(:,1,1)
28659599516SKenneth E. Jansen     &                          + shg(:,2,2) * yl(:,2,1)
28759599516SKenneth E. Jansen     &                          + shg(:,3,2) * yl(:,3,1)
28859599516SKenneth E. Jansen     &                          + shg(:,4,2) * yl(:,4,1)
28959599516SKenneth E. Jansenc
29059599516SKenneth E. Jansen          g2yi(:,2) = g2yi(:,2) + shg(:,1,2) * yl(:,1,2)
29159599516SKenneth E. Jansen     &                          + shg(:,2,2) * yl(:,2,2)
29259599516SKenneth E. Jansen     &                          + shg(:,3,2) * yl(:,3,2)
29359599516SKenneth E. Jansen     &                          + shg(:,4,2) * yl(:,4,2)
29459599516SKenneth E. Jansenc
29559599516SKenneth E. Jansen          g2yi(:,3) = g2yi(:,3) + shg(:,1,2) * yl(:,1,3)
29659599516SKenneth E. Jansen     &                          + shg(:,2,2) * yl(:,2,3)
29759599516SKenneth E. Jansen     &                          + shg(:,3,2) * yl(:,3,3)
29859599516SKenneth E. Jansen     &                          + shg(:,4,2) * yl(:,4,3)
29959599516SKenneth E. Jansenc
30059599516SKenneth E. Jansen          g2yi(:,4) = g2yi(:,4) + shg(:,1,2) * yl(:,1,4)
30159599516SKenneth E. Jansen     &                          + shg(:,2,2) * yl(:,2,4)
30259599516SKenneth E. Jansen     &                          + shg(:,3,2) * yl(:,3,4)
30359599516SKenneth E. Jansen     &                          + shg(:,4,2) * yl(:,4,4)
30459599516SKenneth E. Jansenc
30559599516SKenneth E. Jansen          g2yi(:,5) = g2yi(:,5) + shg(:,1,2) * yl(:,1,5)
30659599516SKenneth E. Jansen     &                          + shg(:,2,2) * yl(:,2,5)
30759599516SKenneth E. Jansen     &                          + shg(:,3,2) * yl(:,3,5)
30859599516SKenneth E. Jansen     &                          + shg(:,4,2) * yl(:,4,5)
30959599516SKenneth E. Jansenc
31059599516SKenneth E. Jansen          g3yi(:,1) = g3yi(:,1) + shg(:,1,3) * yl(:,1,1)
31159599516SKenneth E. Jansen     &                          + shg(:,2,3) * yl(:,2,1)
31259599516SKenneth E. Jansen     &                          + shg(:,3,3) * yl(:,3,1)
31359599516SKenneth E. Jansen     &                          + shg(:,4,3) * yl(:,4,1)
31459599516SKenneth E. Jansenc
31559599516SKenneth E. Jansen          g3yi(:,2) = g3yi(:,2) + shg(:,1,3) * yl(:,1,2)
31659599516SKenneth E. Jansen     &                          + shg(:,2,3) * yl(:,2,2)
31759599516SKenneth E. Jansen     &                          + shg(:,3,3) * yl(:,3,2)
31859599516SKenneth E. Jansen     &                          + shg(:,4,3) * yl(:,4,2)
31959599516SKenneth E. Jansenc
32059599516SKenneth E. Jansen          g3yi(:,3) = g3yi(:,3) + shg(:,1,3) * yl(:,1,3)
32159599516SKenneth E. Jansen     &                          + shg(:,2,3) * yl(:,2,3)
32259599516SKenneth E. Jansen     &                          + shg(:,3,3) * yl(:,3,3)
32359599516SKenneth E. Jansen     &                          + shg(:,4,3) * yl(:,4,3)
32459599516SKenneth E. Jansenc
32559599516SKenneth E. Jansen          g3yi(:,4) = g3yi(:,4) + shg(:,1,3) * yl(:,1,4)
32659599516SKenneth E. Jansen     &                          + shg(:,2,3) * yl(:,2,4)
32759599516SKenneth E. Jansen     &                          + shg(:,3,3) * yl(:,3,4)
32859599516SKenneth E. Jansen     &                          + shg(:,4,3) * yl(:,4,4)
32959599516SKenneth E. Jansenc
33059599516SKenneth E. Jansen          g3yi(:,5) = g3yi(:,5) + shg(:,1,3) * yl(:,1,5)
33159599516SKenneth E. Jansen     &                          + shg(:,2,3) * yl(:,2,5)
33259599516SKenneth E. Jansen     &                          + shg(:,3,3) * yl(:,3,5)
33359599516SKenneth E. Jansen     &                          + shg(:,4,3) * yl(:,4,5)
33459599516SKenneth E. Jansenc
33559599516SKenneth E. Jansen        else
33659599516SKenneth E. Jansen        do n = 1, nshl
33759599516SKenneth E. Jansen          g1yi(:,1) = g1yi(:,1) + shg(:,n,1) * yl(:,n,1)
33859599516SKenneth E. Jansen          g1yi(:,2) = g1yi(:,2) + shg(:,n,1) * yl(:,n,2)
33959599516SKenneth E. Jansen          g1yi(:,3) = g1yi(:,3) + shg(:,n,1) * yl(:,n,3)
34059599516SKenneth E. Jansen          g1yi(:,4) = g1yi(:,4) + shg(:,n,1) * yl(:,n,4)
34159599516SKenneth E. Jansen          g1yi(:,5) = g1yi(:,5) + shg(:,n,1) * yl(:,n,5)
34259599516SKenneth E. Jansenc
34359599516SKenneth E. Jansen          g2yi(:,1) = g2yi(:,1) + shg(:,n,2) * yl(:,n,1)
34459599516SKenneth E. Jansen          g2yi(:,2) = g2yi(:,2) + shg(:,n,2) * yl(:,n,2)
34559599516SKenneth E. Jansen          g2yi(:,3) = g2yi(:,3) + shg(:,n,2) * yl(:,n,3)
34659599516SKenneth E. Jansen          g2yi(:,4) = g2yi(:,4) + shg(:,n,2) * yl(:,n,4)
34759599516SKenneth E. Jansen          g2yi(:,5) = g2yi(:,5) + shg(:,n,2) * yl(:,n,5)
34859599516SKenneth E. Jansenc
34959599516SKenneth E. Jansen          g3yi(:,1) = g3yi(:,1) + shg(:,n,3) * yl(:,n,1)
35059599516SKenneth E. Jansen          g3yi(:,2) = g3yi(:,2) + shg(:,n,3) * yl(:,n,2)
35159599516SKenneth E. Jansen          g3yi(:,3) = g3yi(:,3) + shg(:,n,3) * yl(:,n,3)
35259599516SKenneth E. Jansen          g3yi(:,4) = g3yi(:,4) + shg(:,n,3) * yl(:,n,4)
35359599516SKenneth E. Jansen          g3yi(:,5) = g3yi(:,5) + shg(:,n,3) * yl(:,n,5)
35459599516SKenneth E. Jansenc
35559599516SKenneth E. Jansen        enddo
35659599516SKenneth E. Jansen      endif
35759599516SKenneth E. Jansen
35859599516SKenneth E. Jansen
35959599516SKenneth E. Jansenc
36059599516SKenneth E. Jansenc
36159599516SKenneth E. Jansen         divqi = zero
36259599516SKenneth E. Jansen         idflow = 0
36359599516SKenneth E. Jansen      if (((idiff >= 1) .or. isurf==1).and.
36459599516SKenneth E. Jansen     &     (ires == 3 .or. ires==1)) then
36559599516SKenneth E. Jansen         ttim(32) = ttim(32) - tmr()
36659599516SKenneth E. Jansenc
36759599516SKenneth E. Jansenc     CLASS please ignore all terms related to qi until after you
36859599516SKenneth E. Jansenc     understand EVERYTHING ELSE IN THE CODE.  You may run with
36959599516SKenneth E. Jansenc     idiff=0 for the whole class and everything will be ok never
37059599516SKenneth E. Jansenc     having understood this part.  (leave it for later).
37159599516SKenneth E. Jansenc
37259599516SKenneth E. Jansenc.... compute divergence of diffusive flux vector, qi,i
37359599516SKenneth E. Jansenc
37459599516SKenneth E. Jansen         if(idiff >= 1) then
37559599516SKenneth E. Jansen            idflow = idflow + 4
37659599516SKenneth E. Jansen            do n=1, nshl
37759599516SKenneth E. Jansen
37859599516SKenneth E. Jansen               divqi(:,1) = divqi(:,1) + shg(:,n,1)*ql(:,n,1 )
37959599516SKenneth E. Jansen     &              + shg(:,n,2)*ql(:,n,5 )
38059599516SKenneth E. Jansen     &              + shg(:,n,3)*ql(:,n,9 )
38159599516SKenneth E. Jansen
38259599516SKenneth E. Jansen               divqi(:,2) = divqi(:,2) + shg(:,n,1)*ql(:,n,2 )
38359599516SKenneth E. Jansen     &              + shg(:,n,2)*ql(:,n,6 )
38459599516SKenneth E. Jansen     &              + shg(:,n,3)*ql(:,n,10)
38559599516SKenneth E. Jansen
38659599516SKenneth E. Jansen               divqi(:,3) = divqi(:,3) + shg(:,n,1)*ql(:,n,3 )
38759599516SKenneth E. Jansen     &              + shg(:,n,2)*ql(:,n,7 )
38859599516SKenneth E. Jansen     &              + shg(:,n,3)*ql(:,n,11)
38959599516SKenneth E. Jansen
39059599516SKenneth E. Jansen               divqi(:,4) = divqi(:,4) + shg(:,n,1)*ql(:,n,4 )
39159599516SKenneth E. Jansen     &              + shg(:,n,2)*ql(:,n,8 )
39259599516SKenneth E. Jansen     &              + shg(:,n,3)*ql(:,n,12)
39359599516SKenneth E. Jansen
39459599516SKenneth E. Jansen            enddo
39559599516SKenneth E. Jansen         endif                  !end of idiff
39659599516SKenneth E. Jansenc
39759599516SKenneth E. Jansen         if (isurf .eq. 1) then
39859599516SKenneth E. Jansenc .... divergence of normal calculation
39959599516SKenneth E. Jansen          do n=1, nshl
40059599516SKenneth E. Jansen             divqi(:,idflow+1) = divqi(:,idflow+1)
40159599516SKenneth E. Jansen     &            + shg(:,n,1)*ql(:,n,idflx-2)
40259599516SKenneth E. Jansen     &            + shg(:,n,2)*ql(:,n,idflx-1)
40359599516SKenneth E. Jansen     &            + shg(:,n,3)*ql(:,n,idflx)
40459599516SKenneth E. Jansen          enddo
40559599516SKenneth E. Jansenc .... initialization of some variables
40659599516SKenneth E. Jansen          Sclr = zero
40759599516SKenneth E. Jansen          gradh= zero
40859599516SKenneth E. Jansen          gyti = zero
40959599516SKenneth E. Jansen          sforce=zero
41059599516SKenneth E. Jansen          do i = 1, npro
41159599516SKenneth E. Jansen             do n = 1, nshl
41259599516SKenneth E. Jansen                Sclr(i) = Sclr(i) + shape(i,n) * ycl(i,n,6) !scalar
41359599516SKenneth E. Jansenc
41459599516SKenneth E. Jansenc .... compute the global gradient of Scalar variable
41559599516SKenneth E. Jansenc
41659599516SKenneth E. Jansen                gyti(i,1) = gyti(i,1) + shg(i,n,1) * ycl(i,n,6)
41759599516SKenneth E. Jansen                gyti(i,2) = gyti(i,2) + shg(i,n,2) * ycl(i,n,6)
41859599516SKenneth E. Jansen                gyti(i,3) = gyti(i,3) + shg(i,n,3) * ycl(i,n,6)
41959599516SKenneth E. Jansenc
42059599516SKenneth E. Jansen             enddo
42159599516SKenneth E. Jansen
42259599516SKenneth E. Jansen             if (abs (sclr(i)) .le. epsilon_ls) then
42359599516SKenneth E. Jansen                gradh(i,1) = 0.5/epsilon_ls * (1
42459599516SKenneth E. Jansen     &               + cos(pi*Sclr(i)/epsilon_ls)) * gyti(i,1)
42559599516SKenneth E. Jansen                gradh(i,2) = 0.5/epsilon_ls * (1
42659599516SKenneth E. Jansen     &               + cos(pi*Sclr(i)/epsilon_ls)) * gyti(i,2)
42759599516SKenneth E. Jansen                gradh(i,3) = 0.5/epsilon_ls * (1
42859599516SKenneth E. Jansen     &               + cos(pi*Sclr(i)/epsilon_ls)) * gyti(i,3)
42959599516SKenneth E. Jansen             endif
43059599516SKenneth E. Jansen          enddo                 !end of the loop over npro
43159599516SKenneth E. Jansenc
43259599516SKenneth E. Jansenc .... surface tension force calculation
43359599516SKenneth E. Jansenc .. divide by density now as it gets multiplied in e3res.f, as surface
43459599516SKenneth E. Jansenc    tension force is already in the form of force per unit volume
43559599516SKenneth E. Jansenc
43659599516SKenneth E. Jansen
43759599516SKenneth E. Jansen          weber(:)= Bo
43859599516SKenneth E. Jansen          sforce(:,1) = -(1.0/weber(:)) * divqi(:,idflow+1) !x-direction
43959599516SKenneth E. Jansen     &         *gradh(:,1)/rho(:)
44059599516SKenneth E. Jansen          sforce(:,2) = -(1.0/weber(:)) * divqi(:,idflow+1) !y-direction
44159599516SKenneth E. Jansen     &         *gradh(:,2)/rho(:)
44259599516SKenneth E. Jansen          sforce(:,3) = -(1.0/weber(:)) * divqi(:,idflow+1) !z-direction
44359599516SKenneth E. Jansen     &         *gradh(:,3)/rho(:)
44459599516SKenneth E. Jansen
44559599516SKenneth E. Jansen       endif            ! end of the surface tension force calculation
44659599516SKenneth E. Jansen
44759599516SKenneth E. Jansen         ttim(32) = ttim(32) + secs(0.0)
44859599516SKenneth E. Jansen
44959599516SKenneth E. Jansen      endif                     ! diffusive flux computation
45059599516SKenneth E. Jansenc
45159599516SKenneth E. Jansenc.... return
45259599516SKenneth E. Jansenc
45359599516SKenneth E. Jansen       ttim(20) = ttim(20) + secs(0.0)
45459599516SKenneth E. Jansenc
45559599516SKenneth E. Jansenc.... ----------------------->  dist. kin energy at int. point  <--------------
45659599516SKenneth E. Jansenc
45759599516SKenneth E. Jansenc
45859599516SKenneth E. Jansenc       if (ires .ne. 2 .and. iter.eq.1)  then  !only do at beginning of step
45959599516SKenneth E. Jansenc
46059599516SKenneth E. Jansenc calc exact velocity for a channel at quadrature points.
46159599516SKenneth E. Jansenc
46259599516SKenneth E. Jansenc       dkei=0.0
46359599516SKenneth E. Jansenc
46459599516SKenneth E. Jansenc first guess would be this but it is wrong since it measures the
46559599516SKenneth E. Jansenc error outside of FEM space as well.  This would be correct if we
46659599516SKenneth E. Jansenc wanted to measure error but for disturbance we would like to get
46759599516SKenneth E. Jansenc zero if the solution was nodally exact (i.e no disturbance added to
46859599516SKenneth E. Jansenc the base flow which is nodally exact).
46959599516SKenneth E. Jansenc
47059599516SKenneth E. Jansenc      do n = 1, nshl
47159599516SKenneth E. Jansenc         dkei = dkei + shape(:,n) * xl(:,n,2)  ! this is just the y coord
47259599516SKenneth E. Jansenc      enddo
47359599516SKenneth E. Jansenc         dkei = (1.0-dkei*dkei)  ! u_exact with u_cl=1
47459599516SKenneth E. Jansenc
47559599516SKenneth E. Jansenc
47659599516SKenneth E. Jansenc  correct way
47759599516SKenneth E. Jansenc
47859599516SKenneth E. Jansenc       do n = 1, nshl
47959599516SKenneth E. Jansenc          dkei = dkei + shape(:,n) * (1.0-xl(:,n,2)**2) !u_ex^~  (in FEM space)
48059599516SKenneth E. Jansenc       enddo
48159599516SKenneth E. Jansenc          dkei = (u1-dkei)**2 +u2**2  ! u'^2+v'^2
48259599516SKenneth E. Jansenc          dkei = dkei*WdetJ  ! mult function*W*det of jacobian to
48359599516SKenneth E. Jansenc                              get this quadrature point contribution
48459599516SKenneth E. Jansenc          dke  = dke+sum(dkei) ! we move the sum over elements inside of the
48559599516SKenneth E. Jansenc                              sum over quadrature to save memory (we want
48659599516SKenneth E. Jansenc                              a scalar only)
48759599516SKenneth E. Jansenc       endif
48859599516SKenneth E. Jansen       return
48959599516SKenneth E. Jansen       end
49059599516SKenneth E. Jansen        subroutine e3ivarSclr (ycl,       acl,      acti,
49159599516SKenneth E. Jansen     &                         shape,    shgl,     xl,
49259599516SKenneth E. Jansen     &                         T,        cp,
49359599516SKenneth E. Jansen     &                         dxidx,    Sclr,
494513954efSKenneth E. Jansen     &                         WdetJ,    vort,     gVnrm,
49559599516SKenneth E. Jansen     &                         g1yti,    g2yti,    g3yti,
49659599516SKenneth E. Jansen     &                         rho,      rmu,      con,
49759599516SKenneth E. Jansen     &                         rk,       u1,       u2,
49859599516SKenneth E. Jansen     &                         u3,       shg,      dwl,
49959599516SKenneth E. Jansen     &                         dist2w)
50059599516SKenneth E. Jansenc
50159599516SKenneth E. Jansenc----------------------------------------------------------------------
50259599516SKenneth E. Jansenc
50359599516SKenneth E. Jansenc  This routine computes the variables at integration point.
50459599516SKenneth E. Jansenc
50559599516SKenneth E. Jansenc input:
50659599516SKenneth E. Jansenc  ycl     (npro,nshl,ndof)      : primitive variables
50759599516SKenneth E. Jansenc  actl   (npro,nshl)           : time-deriv of ytl
50859599516SKenneth E. Jansenc  dwl    (npro,nshl)           : distances to wall
50959599516SKenneth E. Jansenc  shape  (npro,nshl)           : element shape-functions
51059599516SKenneth E. Jansenc  shgl   (npro,nsd,nshl)       : element local-grad-shape-functions
51159599516SKenneth E. Jansenc  xl     (npro,nenl,nsd)       : nodal coordinates at current step
51259599516SKenneth E. Jansenc
51359599516SKenneth E. Jansenc output:
51459599516SKenneth E. Jansenc  acti   (npro)                : time-deriv of Scalar variable
51559599516SKenneth E. Jansenc  Sclr   (npro)                : Scalar variable at intpt
51659599516SKenneth E. Jansenc  dist2w (npro)                : distance to nearest wall at intpt
51759599516SKenneth E. Jansenc  WdetJ  (npro)                : weighted Jacobian at intpt
51859599516SKenneth E. Jansenc  vort   (npro)                : vorticity at intpt
519513954efSKenneth E. Jansenc  gVnrm  (npro)                : gradV norm at intpt
52059599516SKenneth E. Jansenc  g1yti  (npro)                : grad-Sclr in direction 1 at intpt
52159599516SKenneth E. Jansenc  g2yti  (npro)                : grad-Sclr in direction 2 at intpt
52259599516SKenneth E. Jansenc  g3yti  (npro)                : grad-Sclr in direction 3 at intpt
52359599516SKenneth E. Jansenc  rho    (npro)                : density at intpt
52459599516SKenneth E. Jansenc  rmu    (npro)                : molecular viscosity
52559599516SKenneth E. Jansenc  con    (npro)                : conductivity
52659599516SKenneth E. Jansenc  rk     (npro)                : kinetic energy
52759599516SKenneth E. Jansenc  u1     (npro)                : x1-velocity component at intpt
52859599516SKenneth E. Jansenc  u2     (npro)                : x2-velocity component at intpt
52959599516SKenneth E. Jansenc  u3     (npro)                : x3-velocity component at intpt
53059599516SKenneth E. Jansenc  shg    (npro,nshl,nsd)       : element global grad-shape-functions at intpt
53159599516SKenneth E. Jansenc
53259599516SKenneth E. Jansenc also used:
53359599516SKenneth E. Jansenc  pres   (npro)                : pressure at intpt
53459599516SKenneth E. Jansenc  T      (npro)                : temperature at intpt
53559599516SKenneth E. Jansenc  ei     (npro)                : internal energy at intpt
53659599516SKenneth E. Jansenc  h      (npro)                : enthalpy at intpt
53759599516SKenneth E. Jansenc  alfap  (npro)                : expansivity at intpt
53859599516SKenneth E. Jansenc  betaT  (npro)                : isothermal compressibility at intpt
53959599516SKenneth E. Jansenc  cp     (npro)                : specific heat at constant pressure at intpt
54059599516SKenneth E. Jansenc  dxidx  (npro,nsd,nsd)        : inverse of deformation gradient at intpt
54159599516SKenneth E. Jansenc  divqi  (npro,nflow-1)         : divergence of diffusive flux
54259599516SKenneth E. Jansenc
54359599516SKenneth E. Jansenc
54459599516SKenneth E. Jansenc Zdenek Johan, Summer 1990. (Modified from e2ivar.f)
54559599516SKenneth E. Jansenc Zdenek Johan, Winter 1991. (Fortran 90)
54659599516SKenneth E. Jansenc Kenneth Jansen, Winter 1997. Primitive Variables
54759599516SKenneth E. Jansenc----------------------------------------------------------------------
54859599516SKenneth E. Jansenc
54959599516SKenneth E. Jansen        include "common.h"
55059599516SKenneth E. Jansenc
55159599516SKenneth E. Jansenc  passed arrays
55259599516SKenneth E. Jansenc
55359599516SKenneth E. Jansen        dimension ycl(npro,nshl,ndof),
55459599516SKenneth E. Jansen     &            acl(npro,nshl,ndof),       acti(npro),
55559599516SKenneth E. Jansen     &            shape(npro,nshl),
55659599516SKenneth E. Jansen     &            shgl(npro,nsd,nshl),       xl(npro,nenl,nsd),
55759599516SKenneth E. Jansen     &            dui(npro,nflow),
55859599516SKenneth E. Jansen     &            g1yi(npro,nflow),
55959599516SKenneth E. Jansen     &            g2yi(npro,nflow),           g3yi(npro,nflow),
56059599516SKenneth E. Jansen     &            shg(npro,nshl,nsd),        dxidx(npro,nsd,nsd),
56159599516SKenneth E. Jansen     &            WdetJ(npro),
56259599516SKenneth E. Jansen     &            rho(npro),                 pres(npro),
56359599516SKenneth E. Jansen     &            T(npro),                   ei(npro),
56459599516SKenneth E. Jansen     &            h(npro),                   alfap(npro),
56559599516SKenneth E. Jansen     &            betaT(npro),               cp(npro),
56659599516SKenneth E. Jansen     &            rk(npro),
56759599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
56859599516SKenneth E. Jansen     &            u3(npro),                  divqi(npro,nflow-1),
56959599516SKenneth E. Jansen     &            ql(npro,nshl,(nflow-1)*nsd),Sclr(npro),
57059599516SKenneth E. Jansen     &            dwl(npro,nenl),
57159599516SKenneth E. Jansen     &            dist2w(npro),
572513954efSKenneth E. Jansen     &            vort(npro),                gVnrm(npro),
57359599516SKenneth E. Jansen     &            rmu(npro),                 con(npro),
57459599516SKenneth E. Jansen     &            g1yti(npro),
57559599516SKenneth E. Jansen     &            g2yti(npro),               g3yti(npro)
57659599516SKenneth E. Jansenc
57759599516SKenneth E. Jansen
57859599516SKenneth E. Jansen        dimension tmp(npro),                  dxdxi(npro,nsd,nsd)
57959599516SKenneth E. Jansenc
58059599516SKenneth E. Jansen        ttim(20) = ttim(20) - tmr()
58159599516SKenneth E. Jansenc
58259599516SKenneth E. Jansenc.... ------------->  Primitive variables at int. point  <--------------
58359599516SKenneth E. Jansenc
58459599516SKenneth E. Jansenc.... compute primitive variables
58559599516SKenneth E. Jansenc
58659599516SKenneth E. Jansen       ttim(11) = ttim(11) - tmr()
58759599516SKenneth E. Jansen
58859599516SKenneth E. Jansen       pres   = zero
58959599516SKenneth E. Jansen       u1     = zero
59059599516SKenneth E. Jansen       u2     = zero
59159599516SKenneth E. Jansen       u3     = zero
59259599516SKenneth E. Jansen       T      = zero
59359599516SKenneth E. Jansen       Sclr   = zero
59459599516SKenneth E. Jansen       dist2w = zero
59559599516SKenneth E. Jansenc
59659599516SKenneth E. Jansen       id = isclr + 5
59759599516SKenneth E. Jansen       do n = 1, nshl
59859599516SKenneth E. Jansenc
59959599516SKenneth E. Jansenc  y(intp)=SUM_{a=1}^nshl (N_a(intp) Ya)
60059599516SKenneth E. Jansenc
60159599516SKenneth E. Jansen          pres   = pres   + shape(:,n) * ycl( :,n,1)
60259599516SKenneth E. Jansen          u1     = u1     + shape(:,n) * ycl( :,n,2)
60359599516SKenneth E. Jansen          u2     = u2     + shape(:,n) * ycl( :,n,3)
60459599516SKenneth E. Jansen          u3     = u3     + shape(:,n) * ycl( :,n,4)
60559599516SKenneth E. Jansen          T      = T      + shape(:,n) * ycl( :,n,5)
60659599516SKenneth E. Jansen          Sclr   = Sclr   + shape(:,n) * ycl(:,n,id)
60759599516SKenneth E. Jansen          if (iRANS.lt.0) then
60859599516SKenneth E. Jansen             dist2w = dist2w + shape(:,n) * dwl(:,n)
60959599516SKenneth E. Jansen          endif
61059599516SKenneth E. Jansen        enddo
61159599516SKenneth E. Jansenc
61259599516SKenneth E. Jansen       ttim(11) = ttim(11) + tmr()
61359599516SKenneth E. Jansenc
61459599516SKenneth E. Jansenc.... ----------------------->  accel. at intp. point  <----------------------
61559599516SKenneth E. Jansenc
61659599516SKenneth E. Jansen
61759599516SKenneth E. Jansen       if (ires .ne. 2)  then
61859599516SKenneth E. Jansen          ttim(12) = ttim(12) - tmr()
61959599516SKenneth E. Jansenc
62059599516SKenneth E. Jansenc.... compute time-derivative of Scalar variable
62159599516SKenneth E. Jansenc
62259599516SKenneth E. Jansen          acti = zero
62359599516SKenneth E. Jansen          do n = 1, nshl
62459599516SKenneth E. Jansen             acti(:)  = acti(:)  + shape(:,n) * acl(:,n,id)
62559599516SKenneth E. Jansen          enddo
62659599516SKenneth E. Jansenc
627*f4d0b58bSKenneth E. Jansen!      flops = flops + 10*nshl*npro
62859599516SKenneth E. Jansen          ttim(12) = ttim(12) + tmr()
62959599516SKenneth E. Jansen       endif                    !ires .ne. 2
63059599516SKenneth E. Jansenc
63159599516SKenneth E. Jansenc.... ----------------->  Thermodynamic Properties  <-------------------
63259599516SKenneth E. Jansenc
63359599516SKenneth E. Jansenc.... compute the kinetic energy
63459599516SKenneth E. Jansenc
63559599516SKenneth E. Jansen        ttim(27) = ttim(27) - tmr()
63659599516SKenneth E. Jansenc
63759599516SKenneth E. Jansen        rk = pt5 * ( u1**2 + u2**2 + u3**2 )
63859599516SKenneth E. Jansenc
63959599516SKenneth E. Jansenc.... get the thermodynamic and material properties
64059599516SKenneth E. Jansenc
64159599516SKenneth E. Jansen        ithm = 7
64259599516SKenneth E. Jansen        call getthm (pres,            T,               Sclr,
64359599516SKenneth E. Jansen     &               rk,              rho,             tmp,
64459599516SKenneth E. Jansen     &               tmp,             tmp,             tmp,
64559599516SKenneth E. Jansen     &               cp,              tmp,             tmp,
64659599516SKenneth E. Jansen     &               tmp,             tmp)
64759599516SKenneth E. Jansenc
64859599516SKenneth E. Jansen        if (iconvsclr.eq.2) rho=one
64959599516SKenneth E. Jansenc
65059599516SKenneth E. Jansen        call getDiffSclr(T,            cp,          rmu,
65159599516SKenneth E. Jansen     &                   tmp,          tmp,         con, rho, Sclr)
65259599516SKenneth E. Jansen
65359599516SKenneth E. Jansen        ttim(27) = ttim(27) + tmr()
65459599516SKenneth E. Jansen
65559599516SKenneth E. Jansen
65659599516SKenneth E. Jansenc
65759599516SKenneth E. Jansenc.... --------------------->  Element Metrics  <-----------------------
65859599516SKenneth E. Jansenc
65959599516SKenneth E. Jansen      call e3metric( xl,         shgl,        dxidx,
66059599516SKenneth E. Jansen     &               shg,        WdetJ)
66159599516SKenneth E. Jansen
66259599516SKenneth E. Jansen
66359599516SKenneth E. Jansen
66459599516SKenneth E. Jansenc.... --------------------->  Global Gradients  <-----------------------
66559599516SKenneth E. Jansenc
66659599516SKenneth E. Jansen        ttim(13) = ttim(13) - tmr()
66759599516SKenneth E. Jansen
66859599516SKenneth E. Jansen        g1yi = zero
66959599516SKenneth E. Jansen        g2yi = zero
67059599516SKenneth E. Jansen        g3yi = zero
67159599516SKenneth E. Jansenc
67259599516SKenneth E. Jansenc.... compute the global gradient of Y-variables
67359599516SKenneth E. Jansenc
67459599516SKenneth E. Jansenc
67559599516SKenneth E. Jansenc  Y_{,x_i}=SUM_{a=1}^nshl (N_{a,x_i}(intp) Ya)
67659599516SKenneth E. Jansenc
67759599516SKenneth E. Jansen        do n = 1, nshl
67859599516SKenneth E. Jansenc         g1yi(:,1) = g1yi(:,1) + shg(:,n,1) * ycl(:,n,1)
679513954efSKenneth E. Jansen          g1yi(:,2) = g1yi(:,2) + shg(:,n,1) * ycl(:,n,2)
68059599516SKenneth E. Jansen          g1yi(:,3) = g1yi(:,3) + shg(:,n,1) * ycl(:,n,3)
68159599516SKenneth E. Jansen          g1yi(:,4) = g1yi(:,4) + shg(:,n,1) * ycl(:,n,4)
68259599516SKenneth E. Jansenc         g1yi(:,5) = g1yi(:,5) + shg(:,n,1) * ycl(:,n,5)
68359599516SKenneth E. Jansenc
68459599516SKenneth E. Jansenc         g2yi(:,1) = g2yi(:,1) + shg(:,n,2) * ycl(:,n,1)
68559599516SKenneth E. Jansen          g2yi(:,2) = g2yi(:,2) + shg(:,n,2) * ycl(:,n,2)
686513954efSKenneth E. Jansen          g2yi(:,3) = g2yi(:,3) + shg(:,n,2) * ycl(:,n,3)
68759599516SKenneth E. Jansen          g2yi(:,4) = g2yi(:,4) + shg(:,n,2) * ycl(:,n,4)
68859599516SKenneth E. Jansenc         g2yi(:,5) = g2yi(:,5) + shg(:,n,2) * ycl(:,n,5)
68959599516SKenneth E. Jansenc
69059599516SKenneth E. Jansenc         g3yi(:,1) = g3yi(:,1) + shg(:,n,3) * ycl(:,n,1)
69159599516SKenneth E. Jansen          g3yi(:,2) = g3yi(:,2) + shg(:,n,3) * ycl(:,n,2)
69259599516SKenneth E. Jansen          g3yi(:,3) = g3yi(:,3) + shg(:,n,3) * ycl(:,n,3)
693513954efSKenneth E. Jansen          g3yi(:,4) = g3yi(:,4) + shg(:,n,3) * ycl(:,n,4)
69459599516SKenneth E. Jansenc         g3yi(:,5) = g3yi(:,5) + shg(:,n,3) * ycl(:,n,5)
69559599516SKenneth E. Jansenc
69659599516SKenneth E. Jansen       enddo
69759599516SKenneth E. Jansen
69859599516SKenneth E. Jansen
69959599516SKenneth E. Jansen
70059599516SKenneth E. Jansen        g1yti = zero
70159599516SKenneth E. Jansen        g2yti = zero
70259599516SKenneth E. Jansen        g3yti = zero
70359599516SKenneth E. Jansenc
70459599516SKenneth E. Jansenc.... compute the global gradient of Scalar variable
70559599516SKenneth E. Jansenc
70659599516SKenneth E. Jansenc
70759599516SKenneth E. Jansenc  Y_{,x_i}=SUM_{a=1}^nshl (N_{a,x_i}(intp) Ya)
70859599516SKenneth E. Jansenc
70959599516SKenneth E. Jansen        do n = 1, nshl
71059599516SKenneth E. Jansen           g1yti(:) = g1yti(:) + shg(:,n,1) * ycl(:,n,id)
71159599516SKenneth E. Jansen           g2yti(:) = g2yti(:) + shg(:,n,2) * ycl(:,n,id)
71259599516SKenneth E. Jansen           g3yti(:) = g3yti(:) + shg(:,n,3) * ycl(:,n,id)
71359599516SKenneth E. Jansenc
71459599516SKenneth E. Jansen        enddo
71559599516SKenneth E. Jansenc **********************************************************
71659599516SKenneth E. Jansenc
71759599516SKenneth E. Jansenc.... compute vorticity at this integration point
71859599516SKenneth E. Jansenc
71959599516SKenneth E. Jansen       vort = sqrt(
72059599516SKenneth E. Jansen     &              (g2yi(:,4)-g3yi(:,3))**2
72159599516SKenneth E. Jansen     &             +(g3yi(:,2)-g1yi(:,4))**2
72259599516SKenneth E. Jansen     &             +(g1yi(:,3)-g2yi(:,2))**2
72359599516SKenneth E. Jansen     &            )
724513954efSKenneth E. Jansen       if(iles.lt.0) then
725513954efSKenneth E. Jansen       gVnrm = sqrt(g1yi(:,2)**2+g2yi(:,2)**2+g3yi(:,2)**2
726513954efSKenneth E. Jansen     &             +g1yi(:,3)**2+g2yi(:,3)**2+g3yi(:,3)**2
727513954efSKenneth E. Jansen     &             +g1yi(:,4)**2+g2yi(:,4)**2+g3yi(:,4)**2)
728513954efSKenneth E. Jansen       endif  ! safe to not define since use is only in same condtnl
729513954efSKenneth E. Jansen
73059599516SKenneth E. Jansenc ***********************************************************
73159599516SKenneth E. Jansen
73259599516SKenneth E. Jansen       ttim(7) = ttim(7) + tmr()
73359599516SKenneth E. Jansen
73459599516SKenneth E. Jansenc
73559599516SKenneth E. Jansenc.... return
73659599516SKenneth E. Jansenc
73759599516SKenneth E. Jansen       ttim(20) = ttim(20) + tmr()
73859599516SKenneth E. Jansen       return
73959599516SKenneth E. Jansen       end
74059599516SKenneth E. Jansen
74159599516SKenneth E. Jansen
742