xref: /phasta/phSolver/compressible/e3bvar.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
159599516SKenneth E. Jansen        subroutine e3bvar (yl,      ycl,     BCB,     shpb,    shglb,
259599516SKenneth E. Jansen     &                     xlb,     lnode,   g1yi,    g2yi,
359599516SKenneth E. Jansen     &                     g3yi,    WdetJb,  bnorm,   pres,    T,
459599516SKenneth E. Jansen     &                     u1,      u2,      u3,      rho,     ei,
559599516SKenneth E. Jansen     &                     cp,      rk,
659599516SKenneth E. Jansen     &                     rou,     p,       tau1n,   tau2n,   tau3n,
7513954efSKenneth E. Jansen     &                     heat,    dNadx)
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansenc----------------------------------------------------------------------
1059599516SKenneth E. Jansenc
1159599516SKenneth E. Jansenc   This routine computes the variables at integration points for
1259599516SKenneth E. Jansenc the boundary element routine.
1359599516SKenneth E. Jansenc
1459599516SKenneth E. Jansenc input:
1559599516SKenneth E. Jansenc  yl     (npro,nshl,nflow)   : primitive variables (perturbed, no scalars)
1659599516SKenneth E. Jansenc  ycl    (npro,nshl,ndof)    : primitive variables
1759599516SKenneth E. Jansenc  BCB    (npro,nshlb,ndBCB)  : Boundary Condition values
1859599516SKenneth E. Jansenc  shpb   (npro,nshl)         : boundary element shape-functions
1959599516SKenneth E. Jansenc  shglb  (npro,nsd,nshl)     : boundary element grad-shape-functions
2059599516SKenneth E. Jansenc  xlb    (npro,nenl,nsd)       : nodal coordinates at current step
2159599516SKenneth E. Jansenc  lnode  (nenb)                : local nodes on the boundary
2259599516SKenneth E. Jansenc
2359599516SKenneth E. Jansenc output:
2459599516SKenneth E. Jansenc  g1yi   (npro,nflow)           : grad-v in direction 1
2559599516SKenneth E. Jansenc  g2yi   (npro,nflow)           : grad-v in direction 2
2659599516SKenneth E. Jansenc  g3yi   (npro,nflow)           : grad-v in direction 3
2759599516SKenneth E. Jansenc  WdetJb (npro)                : weighted Jacobian
2859599516SKenneth E. Jansenc  bnorm  (npro,nsd)            : outward normal
2959599516SKenneth E. Jansenc  pres   (npro)                : pressure
3059599516SKenneth E. Jansenc  T      (npro)                : temperature
3159599516SKenneth E. Jansenc  u1     (npro)                : x1-velocity component
3259599516SKenneth E. Jansenc  u2     (npro)                : x2-velocity component
3359599516SKenneth E. Jansenc  u3     (npro)                : x3-velocity component
3459599516SKenneth E. Jansenc  rho    (npro)                : density
3559599516SKenneth E. Jansenc  ei     (npro)                : internal energy
3659599516SKenneth E. Jansenc  cp     (npro)                : specific energy at constant pressure
3759599516SKenneth E. Jansenc  rk     (npro)                : kinetic energy
3859599516SKenneth E. Jansenc  rou    (npro)                : BC mass flux
3959599516SKenneth E. Jansenc  p      (npro)                : BC pressure
4059599516SKenneth E. Jansenc  tau1n  (npro)                : BC viscous flux 1
4159599516SKenneth E. Jansenc  tau2n  (npro)                : BC viscous flux 2
4259599516SKenneth E. Jansenc  tau3n  (npro)                : BC viscous flux 3
4359599516SKenneth E. Jansenc  heat   (npro)                : BC heat flux
44513954efSKenneth E. Jansenc  dNdx   (npro, nsd)           : BC element shape function gradients
4559599516SKenneth E. Jansenc
4659599516SKenneth E. Jansenc
4759599516SKenneth E. Jansenc Zdenek Johan, Summer 1990.  (Modified from e2bvar.f)
4859599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
4959599516SKenneth E. Jansenc----------------------------------------------------------------------
5059599516SKenneth E. Jansenc
5159599516SKenneth E. Jansen        include "common.h"
5259599516SKenneth E. Jansenc
5359599516SKenneth E. Jansen        dimension yl(npro,nshl,nflow),      BCB(npro,nshlb,ndBCB),
5459599516SKenneth E. Jansen     &            ycl(npro,nshl,ndof),
5559599516SKenneth E. Jansen     &            shpb(npro,nshl),
5659599516SKenneth E. Jansen     &            shglb(npro,nsd,nshl),
5759599516SKenneth E. Jansen     &            xlb(npro,nenl,nsd),
5859599516SKenneth E. Jansen     &            lnode(27),               g1yi(npro,nflow),
5959599516SKenneth E. Jansen     &            g2yi(npro,nflow),        g3yi(npro,nflow),
6059599516SKenneth E. Jansen     &            WdetJb(npro),            bnorm(npro,nsd),
6159599516SKenneth E. Jansen     &            pres(npro),              T(npro),
6259599516SKenneth E. Jansen     &            u1(npro),                u2(npro),
6359599516SKenneth E. Jansen     &            u3(npro),                rho(npro),
6459599516SKenneth E. Jansen     &            ei(npro),                cp(npro),
6559599516SKenneth E. Jansen     &            rk(npro),
6659599516SKenneth E. Jansen     &            rou(npro),               p(npro),
6759599516SKenneth E. Jansen     &            tau1n(npro),             tau2n(npro),
6859599516SKenneth E. Jansen     &            tau3n(npro),             heat(npro)
69513954efSKenneth E. Jansen
7059599516SKenneth E. Jansen        dimension gl1yi(npro,nflow),       gl2yi(npro,nflow),
7159599516SKenneth E. Jansen     &            gl3yi(npro,nflow),       dxdxib(npro,nsd,nsd),
7259599516SKenneth E. Jansen     &            dxidxb(npro,nsd,nsd),    temp(npro),
7359599516SKenneth E. Jansen     &            temp1(npro),             temp2(npro),
74513954efSKenneth E. Jansen     &            temp3(npro),
75513954efSKenneth E. Jansen     &            dNadx(npro, nshl, nsd),  dNadxi(npro, nshl, nsd)
76513954efSKenneth E. Jansen
7759599516SKenneth E. Jansen        dimension h(npro),                 cv(npro),
7859599516SKenneth E. Jansen     &            alfap(npro),             betaT(npro),
7959599516SKenneth E. Jansen     &            gamb(npro),              c(npro),
8059599516SKenneth E. Jansen     &            tmp(npro),
8159599516SKenneth E. Jansen     &            v1(npro,nsd),            v2(npro,nsd)
82513954efSKenneth E. Jansen
83513954efSKenneth E. Jansen        integer   aa
8459599516SKenneth E. Jansenc
8559599516SKenneth E. Jansenc.... ------------------->  integration variables  <--------------------
8659599516SKenneth E. Jansenc
8759599516SKenneth E. Jansenc.... compute the primitive variables at the integration point
8859599516SKenneth E. Jansenc
8959599516SKenneth E. Jansen        pres = zero
9059599516SKenneth E. Jansen        u1   = zero
9159599516SKenneth E. Jansen        u2   = zero
9259599516SKenneth E. Jansen        u3   = zero
9359599516SKenneth E. Jansen        T    = zero
9459599516SKenneth E. Jansenc
9559599516SKenneth E. Jansen        do n = 1, nshlb
9659599516SKenneth E. Jansen          nodlcl = lnode(n)
9759599516SKenneth E. Jansenc
9859599516SKenneth E. Jansen          pres = pres + shpb(:,nodlcl) * yl(:,nodlcl,1)
9959599516SKenneth E. Jansen          u1   = u1   + shpb(:,nodlcl) * yl(:,nodlcl,2)
10059599516SKenneth E. Jansen          u2   = u2   + shpb(:,nodlcl) * yl(:,nodlcl,3)
10159599516SKenneth E. Jansen          u3   = u3   + shpb(:,nodlcl) * yl(:,nodlcl,4)
10259599516SKenneth E. Jansen          T    = T    + shpb(:,nodlcl) * yl(:,nodlcl,5)
10359599516SKenneth E. Jansen        enddo
10459599516SKenneth E. Jansenc
10559599516SKenneth E. Jansenc.... calculate the specific kinetic energy
10659599516SKenneth E. Jansenc
10759599516SKenneth E. Jansen        rk = pt5 * ( u1**2 + u2**2  + u3**2 )
10859599516SKenneth E. Jansenc
10959599516SKenneth E. Jansenc.... get the thermodynamic properties
11059599516SKenneth E. Jansenc
11159599516SKenneth E. Jansen        if (iLSet .ne. 0)then
11259599516SKenneth E. Jansen           temp = zero
11359599516SKenneth E. Jansen           isc=abs(iRANS)+6
11459599516SKenneth E. Jansen           do n = 1, nshlb
11559599516SKenneth E. Jansen              temp = temp + shpb(:,n) * ycl(:,n,isc)
11659599516SKenneth E. Jansen           enddo
11759599516SKenneth E. Jansen        endif
11859599516SKenneth E. Jansen
11959599516SKenneth E. Jansen        ithm = 6
12059599516SKenneth E. Jansen        if (Navier .eq. 1) ithm = 7
12159599516SKenneth E. Jansen        call getthm (pres,            T,                  temp,
12259599516SKenneth E. Jansen     &               rk,              rho,                ei,
12359599516SKenneth E. Jansen     &               h,               tmp,                cv,
12459599516SKenneth E. Jansen     &               cp,              alfap,              betaT,
12559599516SKenneth E. Jansen     &               gamb,            c)
12659599516SKenneth E. Jansenc
12759599516SKenneth E. Jansenc.... ---------------------->  Element Metrics  <-----------------------
12859599516SKenneth E. Jansenc
12959599516SKenneth E. Jansenc.... compute the deformation gradient
13059599516SKenneth E. Jansenc
13159599516SKenneth E. Jansen        dxdxib = zero
13259599516SKenneth E. Jansenc
13359599516SKenneth E. Jansen        do n = 1, nenl
13459599516SKenneth E. Jansen           dxdxib(:,1,1) = dxdxib(:,1,1) + xlb(:,n,1) * shglb(:,1,n)
13559599516SKenneth E. Jansen           dxdxib(:,1,2) = dxdxib(:,1,2) + xlb(:,n,1) * shglb(:,2,n)
13659599516SKenneth E. Jansen           dxdxib(:,1,3) = dxdxib(:,1,3) + xlb(:,n,1) * shglb(:,3,n)
13759599516SKenneth E. Jansen           dxdxib(:,2,1) = dxdxib(:,2,1) + xlb(:,n,2) * shglb(:,1,n)
13859599516SKenneth E. Jansen           dxdxib(:,2,2) = dxdxib(:,2,2) + xlb(:,n,2) * shglb(:,2,n)
13959599516SKenneth E. Jansen           dxdxib(:,2,3) = dxdxib(:,2,3) + xlb(:,n,2) * shglb(:,3,n)
14059599516SKenneth E. Jansen           dxdxib(:,3,1) = dxdxib(:,3,1) + xlb(:,n,3) * shglb(:,1,n)
14159599516SKenneth E. Jansen           dxdxib(:,3,2) = dxdxib(:,3,2) + xlb(:,n,3) * shglb(:,2,n)
14259599516SKenneth E. Jansen           dxdxib(:,3,3) = dxdxib(:,3,3) + xlb(:,n,3) * shglb(:,3,n)
14359599516SKenneth E. Jansen        enddo
14459599516SKenneth E. Jansenc
14559599516SKenneth E. Jansenc.... compute the normal to the boundary
14659599516SKenneth E. Jansenc
14759599516SKenneth E. Jansenc.... compute the normal to the boundary. This is achieved by taking
14859599516SKenneth E. Jansenc     the cross product of two vectors in the plane of the 2-d
14959599516SKenneth E. Jansenc     boundary face.
15059599516SKenneth E. Jansen        v1 = xlb(:,2,:) - xlb(:,1,:)
15159599516SKenneth E. Jansen        v2 = xlb(:,3,:) - xlb(:,1,:)
15259599516SKenneth E. Jansenc
15359599516SKenneth E. Jansenc.....The following are done in order to correct temp1..3
15459599516SKenneth E. Jansenc     based on the results from compressible code.  This is done only
15559599516SKenneth E. Jansenc     for wedges, depending on the bounary face.(tri or quad)
15659599516SKenneth E. Jansenc
15759599516SKenneth E. Jansen        if (lcsyst .eq. 4) then
15859599516SKenneth E. Jansen           temp1 = dxdxib(:,2,1) * dxdxib(:,3,3) -
15959599516SKenneth E. Jansen     &             dxdxib(:,2,3) * dxdxib(:,3,1)
16059599516SKenneth E. Jansen           temp2 = dxdxib(:,3,1) * dxdxib(:,1,3) -
16159599516SKenneth E. Jansen     &             dxdxib(:,3,3) * dxdxib(:,1,1)
16259599516SKenneth E. Jansen           temp3 = dxdxib(:,1,1) * dxdxib(:,2,3) -
16359599516SKenneth E. Jansen     &             dxdxib(:,1,3) * dxdxib(:,2,1)
16459599516SKenneth E. Jansen
16559599516SKenneth E. Jansen        elseif (lcsyst .eq. 1) then
16659599516SKenneth E. Jansen           temp1 = v1(:,2) * v2(:,3) - v2(:,2) * v1(:,3)
16759599516SKenneth E. Jansen           temp2 = v2(:,1) * v1(:,3) - v1(:,1) * v2(:,3)
16859599516SKenneth E. Jansen           temp3 = v1(:,1) * v2(:,2) - v2(:,1) * v1(:,2)
16959599516SKenneth E. Jansen        else
17059599516SKenneth E. Jansen           temp1 = - v1(:,2) * v2(:,3) + v2(:,2) * v1(:,3)
17159599516SKenneth E. Jansen           temp2 = - v2(:,1) * v1(:,3) + v1(:,1) * v2(:,3)
17259599516SKenneth E. Jansen           temp3 = - v1(:,1) * v2(:,2) + v2(:,1) * v1(:,2)
17359599516SKenneth E. Jansen        endif
17459599516SKenneth E. Jansenc
17559599516SKenneth E. Jansen        temp       = one / sqrt ( temp1**2 + temp2**2 + temp3**2 )
17659599516SKenneth E. Jansen        bnorm(:,1) = temp1 * temp
17759599516SKenneth E. Jansen        bnorm(:,2) = temp2 * temp
17859599516SKenneth E. Jansen        bnorm(:,3) = temp3 * temp
17959599516SKenneth E. Jansenc
18059599516SKenneth E. Jansen
18159599516SKenneth E. Jansen        if (lcsyst .eq. 3) then
18259599516SKenneth E. Jansen           WdetJb     = (1 - Qwtb(lcsyst,intp)) / (four*temp)
18359599516SKenneth E. Jansen        elseif (lcsyst .eq. 4) then
18459599516SKenneth E. Jansenc
18559599516SKenneth E. Jansenc  quad face wedges have a conflict in lnode ordering that makes the
18659599516SKenneth E. Jansenc  normal negative
18759599516SKenneth E. Jansenc
18859599516SKenneth E. Jansenc          bnorm=-bnorm
18959599516SKenneth E. Jansenc
19059599516SKenneth E. Jansen          WdetJb     = Qwtb(lcsyst,intp) / temp
19159599516SKenneth E. Jansen        else
19259599516SKenneth E. Jansen           WdetJb     = Qwtb(lcsyst,intp) / (four*temp)
19359599516SKenneth E. Jansen        endif
19459599516SKenneth E. Jansenc
19559599516SKenneth E. Jansenc.... -------------------------->  Grad-V  <----------------------------
19659599516SKenneth E. Jansenc
19759599516SKenneth E. Jansenc.... compute grad-v for Navier-Stokes terms
19859599516SKenneth E. Jansenc
19959599516SKenneth E. Jansen        if (Navier .eq. 1) then
20059599516SKenneth E. Jansenc
20159599516SKenneth E. Jansenc.... compute the inverse of deformation gradient
20259599516SKenneth E. Jansenc
20359599516SKenneth E. Jansen          dxidxb(:,1,1) =   dxdxib(:,2,2) * dxdxib(:,3,3)
20459599516SKenneth E. Jansen     &                    - dxdxib(:,3,2) * dxdxib(:,2,3)
20559599516SKenneth E. Jansen          dxidxb(:,1,2) =   dxdxib(:,3,2) * dxdxib(:,1,3)
20659599516SKenneth E. Jansen     &                    - dxdxib(:,1,2) * dxdxib(:,3,3)
20759599516SKenneth E. Jansen          dxidxb(:,1,3) =   dxdxib(:,1,2) * dxdxib(:,2,3)
20859599516SKenneth E. Jansen     &                    - dxdxib(:,1,3) * dxdxib(:,2,2)
20959599516SKenneth E. Jansen          temp          = one / ( dxidxb(:,1,1) * dxdxib(:,1,1)
21059599516SKenneth E. Jansen     &                          + dxidxb(:,1,2) * dxdxib(:,2,1)
21159599516SKenneth E. Jansen     &                          + dxidxb(:,1,3) * dxdxib(:,3,1) )
21259599516SKenneth E. Jansen          dxidxb(:,1,1) =  dxidxb(:,1,1) * temp
21359599516SKenneth E. Jansen          dxidxb(:,1,2) =  dxidxb(:,1,2) * temp
21459599516SKenneth E. Jansen          dxidxb(:,1,3) =  dxidxb(:,1,3) * temp
21559599516SKenneth E. Jansen          dxidxb(:,2,1) = (dxdxib(:,2,3) * dxdxib(:,3,1)
21659599516SKenneth E. Jansen     &                   - dxdxib(:,2,1) * dxdxib(:,3,3)) * temp
21759599516SKenneth E. Jansen          dxidxb(:,2,2) = (dxdxib(:,1,1) * dxdxib(:,3,3)
21859599516SKenneth E. Jansen     &                   - dxdxib(:,3,1) * dxdxib(:,1,3)) * temp
21959599516SKenneth E. Jansen          dxidxb(:,2,3) = (dxdxib(:,2,1) * dxdxib(:,1,3)
22059599516SKenneth E. Jansen     &                   - dxdxib(:,1,1) * dxdxib(:,2,3)) * temp
22159599516SKenneth E. Jansen          dxidxb(:,3,1) = (dxdxib(:,2,1) * dxdxib(:,3,2)
22259599516SKenneth E. Jansen     &                   - dxdxib(:,2,2) * dxdxib(:,3,1)) * temp
22359599516SKenneth E. Jansen          dxidxb(:,3,2) = (dxdxib(:,3,1) * dxdxib(:,1,2)
22459599516SKenneth E. Jansen     &                   - dxdxib(:,1,1) * dxdxib(:,3,2)) * temp
22559599516SKenneth E. Jansen          dxidxb(:,3,3) = (dxdxib(:,1,1) * dxdxib(:,2,2)
22659599516SKenneth E. Jansen     &                   - dxdxib(:,1,2) * dxdxib(:,2,1)) * temp
22759599516SKenneth E. Jansenc
22859599516SKenneth E. Jansenc.... compute local-grad-Y
22959599516SKenneth E. Jansenc
23059599516SKenneth E. Jansen          gl1yi = zero
23159599516SKenneth E. Jansen          gl2yi = zero
23259599516SKenneth E. Jansen          gl3yi = zero
23359599516SKenneth E. Jansenc
23459599516SKenneth E. Jansen          do n = 1, nshl
23559599516SKenneth E. Jansen            gl1yi(:,1) = gl1yi(:,1) + shglb(:,1,n) * yl(:,n,1)
23659599516SKenneth E. Jansen            gl1yi(:,2) = gl1yi(:,2) + shglb(:,1,n) * yl(:,n,2)
23759599516SKenneth E. Jansen            gl1yi(:,3) = gl1yi(:,3) + shglb(:,1,n) * yl(:,n,3)
23859599516SKenneth E. Jansen            gl1yi(:,4) = gl1yi(:,4) + shglb(:,1,n) * yl(:,n,4)
23959599516SKenneth E. Jansen            gl1yi(:,5) = gl1yi(:,5) + shglb(:,1,n) * yl(:,n,5)
24059599516SKenneth E. Jansenc
24159599516SKenneth E. Jansen            gl2yi(:,1) = gl2yi(:,1) + shglb(:,2,n) * yl(:,n,1)
24259599516SKenneth E. Jansen            gl2yi(:,2) = gl2yi(:,2) + shglb(:,2,n) * yl(:,n,2)
24359599516SKenneth E. Jansen            gl2yi(:,3) = gl2yi(:,3) + shglb(:,2,n) * yl(:,n,3)
24459599516SKenneth E. Jansen            gl2yi(:,4) = gl2yi(:,4) + shglb(:,2,n) * yl(:,n,4)
24559599516SKenneth E. Jansen            gl2yi(:,5) = gl2yi(:,5) + shglb(:,2,n) * yl(:,n,5)
24659599516SKenneth E. Jansenc
24759599516SKenneth E. Jansen            gl3yi(:,1) = gl3yi(:,1) + shglb(:,3,n) * yl(:,n,1)
24859599516SKenneth E. Jansen            gl3yi(:,2) = gl3yi(:,2) + shglb(:,3,n) * yl(:,n,2)
24959599516SKenneth E. Jansen            gl3yi(:,3) = gl3yi(:,3) + shglb(:,3,n) * yl(:,n,3)
25059599516SKenneth E. Jansen            gl3yi(:,4) = gl3yi(:,4) + shglb(:,3,n) * yl(:,n,4)
25159599516SKenneth E. Jansen            gl3yi(:,5) = gl3yi(:,5) + shglb(:,3,n) * yl(:,n,5)
25259599516SKenneth E. Jansen          enddo
25359599516SKenneth E. Jansenc
25459599516SKenneth E. Jansenc.... convert local-grads to global-grads
25559599516SKenneth E. Jansenc
25659599516SKenneth E. Jansen          g1yi(:,2) = dxidxb(:,1,1) * gl1yi(:,2) +
25759599516SKenneth E. Jansen     &                dxidxb(:,2,1) * gl2yi(:,2) +
25859599516SKenneth E. Jansen     &                dxidxb(:,3,1) * gl3yi(:,2)
25959599516SKenneth E. Jansen          g2yi(:,2) = dxidxb(:,1,2) * gl1yi(:,2) +
26059599516SKenneth E. Jansen     &                dxidxb(:,2,2) * gl2yi(:,2) +
26159599516SKenneth E. Jansen     &                dxidxb(:,3,2) * gl3yi(:,2)
26259599516SKenneth E. Jansen          g3yi(:,2) = dxidxb(:,1,3) * gl1yi(:,2) +
26359599516SKenneth E. Jansen     &                dxidxb(:,2,3) * gl2yi(:,2) +
26459599516SKenneth E. Jansen     &                dxidxb(:,3,3) * gl3yi(:,2)
26559599516SKenneth E. Jansenc
26659599516SKenneth E. Jansen          g1yi(:,3) = dxidxb(:,1,1) * gl1yi(:,3) +
26759599516SKenneth E. Jansen     &                dxidxb(:,2,1) * gl2yi(:,3) +
26859599516SKenneth E. Jansen     &                dxidxb(:,3,1) * gl3yi(:,3)
26959599516SKenneth E. Jansen          g2yi(:,3) = dxidxb(:,1,2) * gl1yi(:,3) +
27059599516SKenneth E. Jansen     &                dxidxb(:,2,2) * gl2yi(:,3) +
27159599516SKenneth E. Jansen     &                dxidxb(:,3,2) * gl3yi(:,3)
27259599516SKenneth E. Jansen          g3yi(:,3) = dxidxb(:,1,3) * gl1yi(:,3) +
27359599516SKenneth E. Jansen     &                dxidxb(:,2,3) * gl2yi(:,3) +
27459599516SKenneth E. Jansen     &                dxidxb(:,3,3) * gl3yi(:,3)
27559599516SKenneth E. Jansenc
27659599516SKenneth E. Jansen          g1yi(:,4) = dxidxb(:,1,1) * gl1yi(:,4) +
27759599516SKenneth E. Jansen     &                dxidxb(:,2,1) * gl2yi(:,4) +
27859599516SKenneth E. Jansen     &                dxidxb(:,3,1) * gl3yi(:,4)
27959599516SKenneth E. Jansen          g2yi(:,4) = dxidxb(:,1,2) * gl1yi(:,4) +
28059599516SKenneth E. Jansen     &                dxidxb(:,2,2) * gl2yi(:,4) +
28159599516SKenneth E. Jansen     &                dxidxb(:,3,2) * gl3yi(:,4)
28259599516SKenneth E. Jansen          g3yi(:,4) = dxidxb(:,1,3) * gl1yi(:,4) +
28359599516SKenneth E. Jansen     &                dxidxb(:,2,3) * gl2yi(:,4) +
28459599516SKenneth E. Jansen     &                dxidxb(:,3,3) * gl3yi(:,4)
28559599516SKenneth E. Jansenc
28659599516SKenneth E. Jansen          g1yi(:,5) = dxidxb(:,1,1) * gl1yi(:,5) +
28759599516SKenneth E. Jansen     &                dxidxb(:,2,1) * gl2yi(:,5) +
28859599516SKenneth E. Jansen     &                dxidxb(:,3,1) * gl3yi(:,5)
28959599516SKenneth E. Jansen          g2yi(:,5) = dxidxb(:,1,2) * gl1yi(:,5) +
29059599516SKenneth E. Jansen     &                dxidxb(:,2,2) * gl2yi(:,5) +
29159599516SKenneth E. Jansen     &                dxidxb(:,3,2) * gl3yi(:,5)
29259599516SKenneth E. Jansen          g3yi(:,5) = dxidxb(:,1,3) * gl1yi(:,5) +
29359599516SKenneth E. Jansen     &                dxidxb(:,2,3) * gl2yi(:,5) +
29459599516SKenneth E. Jansen     &                dxidxb(:,3,3) * gl3yi(:,5)
29559599516SKenneth E. Jansenc
29659599516SKenneth E. Jansenc.... end grad-v
29759599516SKenneth E. Jansenc
298513954efSKenneth E. Jansen          !Compute the gradient of the shape function for heat flux's
299513954efSKenneth E. Jansen          !contribution to lhsk
300513954efSKenneth E. Jansen          if(iLHScond > 0) then
301513954efSKenneth E. Jansen            dNadx = zero
302513954efSKenneth E. Jansen
303513954efSKenneth E. Jansen            !dNdx(a,i) = dN_a / dx_i
304513954efSKenneth E. Jansen
305513954efSKenneth E. Jansen            do aa = 1, nshl  !TODO: get rid of the intermediary dNadxi
306513954efSKenneth E. Jansen                             !shglb(:,nsd,a=1)* N(:,a=1)
307513954efSKenneth E. Jansen              dNadxi(:,aa,1) = shglb(:,1,aa) * 1 !would normally be a sum over
308513954efSKenneth E. Jansen              dNadxi(:,aa,2) = shglb(:,2,aa) * 1 !all nodes, but N = 0 for a /= 1
309513954efSKenneth E. Jansen              dNadxi(:,aa,3) = shglb(:,3,aa) * 1
310513954efSKenneth E. Jansen            enddo
311513954efSKenneth E. Jansen
312513954efSKenneth E. Jansen            do aa = 1, nshl
313513954efSKenneth E. Jansen              dNadx(:,aa,1) = dNadxi(:,aa,1) * dxidxb(:,1,1) +
314513954efSKenneth E. Jansen     &                        dNadxi(:,aa,2) * dxidxb(:,2,1) +
315513954efSKenneth E. Jansen     &                        dNadxi(:,aa,3) * dxidxb(:,3,1)
316513954efSKenneth E. Jansen              dNadx(:,aa,2) = dNadxi(:,aa,1) * dxidxb(:,1,2) +
317513954efSKenneth E. Jansen     &                        dNadxi(:,aa,2) * dxidxb(:,2,2) +
318513954efSKenneth E. Jansen     &                        dNadxi(:,aa,3) * dxidxb(:,3,2)
319513954efSKenneth E. Jansen              dNadx(:,aa,3) = dNadxi(:,aa,1) * dxidxb(:,1,3) +
320513954efSKenneth E. Jansen     &                        dNadxi(:,aa,2) * dxidxb(:,2,3) +
321513954efSKenneth E. Jansen     &                        dNadxi(:,aa,3) * dxidxb(:,3,3)
322513954efSKenneth E. Jansen            enddo
323513954efSKenneth E. Jansen          endif
324513954efSKenneth E. Jansen
32559599516SKenneth E. Jansen        endif
32659599516SKenneth E. Jansenc
32759599516SKenneth E. Jansenc.... -------------------->  Boundary Conditions  <--------------------
32859599516SKenneth E. Jansenc
32959599516SKenneth E. Jansenc.... compute the Euler boundary conditions
33059599516SKenneth E. Jansenc
33159599516SKenneth E. Jansen        rou = zero
33259599516SKenneth E. Jansen        p   = zero
33359599516SKenneth E. Jansenc
33459599516SKenneth E. Jansen        do n = 1, nshlb
33559599516SKenneth E. Jansen          nodlcl = lnode(n)
33659599516SKenneth E. Jansenc
33759599516SKenneth E. Jansen          rou = rou + shpb(:,nodlcl) * BCB(:,n,1)
33859599516SKenneth E. Jansen          p   = p   + shpb(:,nodlcl) * BCB(:,n,2)
33959599516SKenneth E. Jansen        enddo
34059599516SKenneth E. Jansenc
34159599516SKenneth E. Jansenc.... compute the Navier-Stokes boundary conditions
34259599516SKenneth E. Jansenc
34359599516SKenneth E. Jansen        if (Navier .eq. 1) then
34459599516SKenneth E. Jansenc
34559599516SKenneth E. Jansen          tau1n = zero
34659599516SKenneth E. Jansen          tau2n = zero
34759599516SKenneth E. Jansen          tau3n = zero
34859599516SKenneth E. Jansen          heat  = zero
34959599516SKenneth E. Jansenc
35059599516SKenneth E. Jansen          do n = 1, nshlb
35159599516SKenneth E. Jansen            nodlcl = lnode(n)
35259599516SKenneth E. Jansenc
35359599516SKenneth E. Jansen            tau1n = tau1n + shpb(:,nodlcl) * BCB(:,n,3)
35459599516SKenneth E. Jansen            tau2n = tau2n + shpb(:,nodlcl) * BCB(:,n,4)
35559599516SKenneth E. Jansen            tau3n = tau3n + shpb(:,nodlcl) * BCB(:,n,5)
35659599516SKenneth E. Jansen            heat  = heat  + shpb(:,nodlcl) * BCB(:,n,6)
35759599516SKenneth E. Jansen          enddo
35859599516SKenneth E. Jansenc
35959599516SKenneth E. Jansenc.... flop count
36059599516SKenneth E. Jansenc
361*f4d0b58bSKenneth E. Jansen!      flops = flops + (184+30*nshl+8*nshlb)*npro
36259599516SKenneth E. Jansenc
36359599516SKenneth E. Jansen        endif
36459599516SKenneth E. Jansenc
36559599516SKenneth E. Jansenc.... flop count
36659599516SKenneth E. Jansenc
367*f4d0b58bSKenneth E. Jansen!      flops = flops + (27+18*nshl+14*nshlb)*npro
36859599516SKenneth E. Jansenc
36959599516SKenneth E. Jansenc.... return
37059599516SKenneth E. Jansenc
37159599516SKenneth E. Jansen        return
37259599516SKenneth E. Jansen        end
37359599516SKenneth E. Jansenc
37459599516SKenneth E. Jansenc
37559599516SKenneth E. Jansenc
37659599516SKenneth E. Jansen        subroutine e3bvarSclr(ycl,      BCB,     shpb,    shglb,
37759599516SKenneth E. Jansen     &                        xlb,     lnode,
37859599516SKenneth E. Jansen     &                        u1,      u2,      u3,
37959599516SKenneth E. Jansen     &                        g1yti,   g2yti,   g3yti,   WdetJb,
38059599516SKenneth E. Jansen     &                        bnorm,   T,       rho,     cp,      rou,
38159599516SKenneth E. Jansen     &                        Sclr,    SclrF)
38259599516SKenneth E. Jansenc
38359599516SKenneth E. Jansenc----------------------------------------------------------------------
38459599516SKenneth E. Jansenc
38559599516SKenneth E. Jansenc   This routine computes the variables at integration points for
38659599516SKenneth E. Jansenc the boundary element routine.
38759599516SKenneth E. Jansenc
38859599516SKenneth E. Jansenc input:
38959599516SKenneth E. Jansenc  ycl     (npro,nshl,ndof)      : Y variables
39059599516SKenneth E. Jansenc  BCB    (npro,nenbl,ndBCB)    : Boundary Condition values
39159599516SKenneth E. Jansenc  shpb   (npro,nen)            : boundary element shape-functions
39259599516SKenneth E. Jansenc  shglb  (nsd,nen)             : boundary element grad-shape-functions
39359599516SKenneth E. Jansenc  xlb    (npro,nshl,nsd)       : nodal coordinates at current step
39459599516SKenneth E. Jansenc  lnode  (nenb)                : local nodes on the boundary
39559599516SKenneth E. Jansenc
39659599516SKenneth E. Jansenc output:
39759599516SKenneth E. Jansenc  g1yti  (npro)
39859599516SKenneth E. Jansenc  g2yti  (npro)
39959599516SKenneth E. Jansenc  g3yti  (npro)
40059599516SKenneth E. Jansenc  WdetJb (npro)                : weighted Jacobian
40159599516SKenneth E. Jansenc  bnorm  (npro,nsd)            : outward normal
40259599516SKenneth E. Jansenc  T      (npro)                : temperature
40359599516SKenneth E. Jansenc  rho    (npro)                : density
40459599516SKenneth E. Jansenc  cp     (npro)                : specific energy at constant pressure
40559599516SKenneth E. Jansenc  rou    (npro)                : BC mass flux
40659599516SKenneth E. Jansenc  SclrF  (npro)                : BC Scalar  flux
40759599516SKenneth E. Jansenc
40859599516SKenneth E. Jansenc Zdenek Johan, Summer 1990.  (Modified from e2bvar.f)
40959599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
41059599516SKenneth E. Jansenc----------------------------------------------------------------------
41159599516SKenneth E. Jansenc
41259599516SKenneth E. Jansen        include "common.h"
41359599516SKenneth E. Jansenc
41459599516SKenneth E. Jansen        dimension ycl(npro,nshl,ndof),        BCB(npro,nshlb,ndBCB),
41559599516SKenneth E. Jansen     &            shpb(npro,nshl),           shglb(npro,nsd,nshl),
41659599516SKenneth E. Jansen     &            xlb(npro,nshl,nsd),
41759599516SKenneth E. Jansen     &            lnode(27),
41859599516SKenneth E. Jansen     &            g1yti(npro),               g2yti(npro),
41959599516SKenneth E. Jansen     &            g3yti(npro),
42059599516SKenneth E. Jansen     &            WdetJb(npro),              bnorm(npro,nsd),
42159599516SKenneth E. Jansen     &            pres(npro),                T(npro),
42259599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
42359599516SKenneth E. Jansen     &            u3(npro),                  rho(npro),
42459599516SKenneth E. Jansen     &            ei(npro),                  cp(npro),
42559599516SKenneth E. Jansen     &            rk(npro),                  Sclr(npro),
42659599516SKenneth E. Jansen     &            rou(npro),
42759599516SKenneth E. Jansen     &            SclrF(npro)
42859599516SKenneth E. Jansenc
42959599516SKenneth E. Jansen        dimension dxdxib(npro,nsd,nsd),
43059599516SKenneth E. Jansen     &            dxidxb(npro,nsd,nsd),      temp(npro),
43159599516SKenneth E. Jansen     &            temp1(npro),               temp2(npro),
43259599516SKenneth E. Jansen     &            temp3(npro),               gl1yti(npro),
43359599516SKenneth E. Jansen     &            gl2yti(npro),              gl3yti(npro)
43459599516SKenneth E. Jansenc
43559599516SKenneth E. Jansen        dimension h(npro),                   cv(npro),
43659599516SKenneth E. Jansen     &            alfap(npro),               betaT(npro),
43759599516SKenneth E. Jansen     &            gamb(npro),                c(npro),
43859599516SKenneth E. Jansen     &            tmp(npro),                 v1(npro,nsd),
43959599516SKenneth E. Jansen     &            v2(npro,nsd)
44059599516SKenneth E. Jansenc
44159599516SKenneth E. Jansenc.... ------------------->  integration variables  <--------------------
44259599516SKenneth E. Jansenc
44359599516SKenneth E. Jansenc.... compute the primitive variables at the integration point
44459599516SKenneth E. Jansenc
44559599516SKenneth E. Jansen        pres = zero
44659599516SKenneth E. Jansen        u1   = zero
44759599516SKenneth E. Jansen        u2   = zero
44859599516SKenneth E. Jansen        u3   = zero
44959599516SKenneth E. Jansen        T    = zero
45059599516SKenneth E. Jansen        Sclr = zero
45159599516SKenneth E. Jansen
45259599516SKenneth E. Jansen        id  = isclr+5
45359599516SKenneth E. Jansen        ibb = isclr+6
45459599516SKenneth E. Jansenc
45559599516SKenneth E. Jansen        do n = 1, nshlb
45659599516SKenneth E. Jansen          nodlcl = lnode(n)
45759599516SKenneth E. Jansenc
45859599516SKenneth E. Jansen          pres = pres + shpb(:,nodlcl) * ycl(:,nodlcl,1)
45959599516SKenneth E. Jansen          u1   = u1   + shpb(:,nodlcl) * ycl(:,nodlcl,2)
46059599516SKenneth E. Jansen          u2   = u2   + shpb(:,nodlcl) * ycl(:,nodlcl,3)
46159599516SKenneth E. Jansen          u3   = u3   + shpb(:,nodlcl) * ycl(:,nodlcl,4)
46259599516SKenneth E. Jansen          T    = T    + shpb(:,nodlcl) * ycl(:,nodlcl,5)
46359599516SKenneth E. Jansen          Sclr = Sclr + shpb(:,nodlcl) * ycl(:,nodlcl,id)
46459599516SKenneth E. Jansen        enddo
46559599516SKenneth E. Jansenc
46659599516SKenneth E. Jansenc.... calculate the specific kinetic energy
46759599516SKenneth E. Jansenc
46859599516SKenneth E. Jansen        rk = pt5 * ( u1**2 + u2**2  + u3**2 )
46959599516SKenneth E. Jansenc
47059599516SKenneth E. Jansenc.... get the thermodynamic properties
47159599516SKenneth E. Jansenc
47259599516SKenneth E. Jansen        ithm = 6
47359599516SKenneth E. Jansen        if (Navier .eq. 1) ithm = 7
47459599516SKenneth E. Jansen        call getthm (pres,            T,                  Sclr,
47559599516SKenneth E. Jansen     &               rk,              rho,                ei,
47659599516SKenneth E. Jansen     &               h,               tmp,                cv,
47759599516SKenneth E. Jansen     &               cp,              alfap,              betaT,
47859599516SKenneth E. Jansen     &               gamb,            c)
47959599516SKenneth E. Jansenc
48059599516SKenneth E. Jansen       if (iconvsclr.eq.2) rho=one
48159599516SKenneth E. Jansenc
48259599516SKenneth E. Jansenc.... ---------------------->  Element Metrics  <-----------------------
48359599516SKenneth E. Jansenc
48459599516SKenneth E. Jansenc.... compute the deformation gradient
48559599516SKenneth E. Jansenc
48659599516SKenneth E. Jansen        dxdxib = zero
48759599516SKenneth E. Jansenc
48859599516SKenneth E. Jansen        do n = 1, nshl
48959599516SKenneth E. Jansen           dxdxib(:,1,1) = dxdxib(:,1,1) + xlb(:,n,1) * shglb(:,1,n)
49059599516SKenneth E. Jansen           dxdxib(:,1,2) = dxdxib(:,1,2) + xlb(:,n,1) * shglb(:,2,n)
49159599516SKenneth E. Jansen           dxdxib(:,1,3) = dxdxib(:,1,3) + xlb(:,n,1) * shglb(:,3,n)
49259599516SKenneth E. Jansen           dxdxib(:,2,1) = dxdxib(:,2,1) + xlb(:,n,2) * shglb(:,1,n)
49359599516SKenneth E. Jansen           dxdxib(:,2,2) = dxdxib(:,2,2) + xlb(:,n,2) * shglb(:,2,n)
49459599516SKenneth E. Jansen           dxdxib(:,2,3) = dxdxib(:,2,3) + xlb(:,n,2) * shglb(:,3,n)
49559599516SKenneth E. Jansen           dxdxib(:,3,1) = dxdxib(:,3,1) + xlb(:,n,3) * shglb(:,1,n)
49659599516SKenneth E. Jansen           dxdxib(:,3,2) = dxdxib(:,3,2) + xlb(:,n,3) * shglb(:,2,n)
49759599516SKenneth E. Jansen           dxdxib(:,3,3) = dxdxib(:,3,3) + xlb(:,n,3) * shglb(:,3,n)
49859599516SKenneth E. Jansen        enddo
49959599516SKenneth E. Jansenc
50059599516SKenneth E. Jansenc.... compute the normal to the boundary
50159599516SKenneth E. Jansenc
50259599516SKenneth E. Jansenc
50359599516SKenneth E. Jansen        v1 = xlb(:,2,:) - xlb(:,1,:)
50459599516SKenneth E. Jansen        v2 = xlb(:,3,:) - xlb(:,1,:)
50559599516SKenneth E. Jansenc
50659599516SKenneth E. Jansenc.... compute the normal to the boundary. This is achieved by taking
50759599516SKenneth E. Jansenc     the cross product of two vectors in the plane of the 2-d
50859599516SKenneth E. Jansenc     boundary face.
50959599516SKenneth E. Jansenc
51059599516SKenneth E. Jansenc.....The following are done in order to correct temp1..3
51159599516SKenneth E. Jansenc     based on the results from compressible code.  This is done only
51259599516SKenneth E. Jansenc     for wedges, depending on the bounary face.(tri or quad)
51359599516SKenneth E. Jansenc
51459599516SKenneth E. Jansen        if (lcsyst .eq. 4) then
51559599516SKenneth E. Jansen           temp1 = dxdxib(:,2,1) * dxdxib(:,3,3) -
51659599516SKenneth E. Jansen     &             dxdxib(:,2,3) * dxdxib(:,3,1)
51759599516SKenneth E. Jansen           temp2 = dxdxib(:,3,1) * dxdxib(:,1,3) -
51859599516SKenneth E. Jansen     &             dxdxib(:,3,3) * dxdxib(:,1,1)
51959599516SKenneth E. Jansen           temp3 = dxdxib(:,1,1) * dxdxib(:,2,3) -
52059599516SKenneth E. Jansen     &             dxdxib(:,1,3) * dxdxib(:,2,1)
52159599516SKenneth E. Jansen
52259599516SKenneth E. Jansen        elseif (lcsyst .eq. 1) then
52359599516SKenneth E. Jansen           temp1 = v1(:,2) * v2(:,3) - v2(:,2) * v1(:,3)
52459599516SKenneth E. Jansen           temp2 = v2(:,1) * v1(:,3) - v1(:,1) * v2(:,3)
52559599516SKenneth E. Jansen           temp3 = v1(:,1) * v2(:,2) - v2(:,1) * v1(:,2)
52659599516SKenneth E. Jansen        else
52759599516SKenneth E. Jansen           temp1 = - v1(:,2) * v2(:,3) + v2(:,2) * v1(:,3)
52859599516SKenneth E. Jansen           temp2 = - v2(:,1) * v1(:,3) + v1(:,1) * v2(:,3)
52959599516SKenneth E. Jansen           temp3 = - v1(:,1) * v2(:,2) + v2(:,1) * v1(:,2)
53059599516SKenneth E. Jansen        endif
53159599516SKenneth E. Jansenc
53259599516SKenneth E. Jansen        temp       = one / sqrt ( temp1**2 + temp2**2 + temp3**2 )
53359599516SKenneth E. Jansen        bnorm(:,1) = temp1 * temp
53459599516SKenneth E. Jansen        bnorm(:,2) = temp2 * temp
53559599516SKenneth E. Jansen        bnorm(:,3) = temp3 * temp
53659599516SKenneth E. Jansenc
53759599516SKenneth E. Jansen
53859599516SKenneth E. Jansen        if (lcsyst .eq. 3) then
53959599516SKenneth E. Jansen           WdetJb     = (1 - Qwtb(lcsyst,intp)) / (four*temp)
54059599516SKenneth E. Jansen        elseif (lcsyst .eq. 4) then
54159599516SKenneth E. Jansen          WdetJb     = Qwtb(lcsyst,intp)/ temp
54259599516SKenneth E. Jansen        else
54359599516SKenneth E. Jansen           WdetJb     =Qwtb(lcsyst,intp) / (four*temp)
54459599516SKenneth E. Jansen        endif
54559599516SKenneth E. Jansenc
54659599516SKenneth E. Jansenc.... -------------------------->  Grad-V  <----------------------------
54759599516SKenneth E. Jansenc
54859599516SKenneth E. Jansenc.... compute grad-v for Navier-Stokes terms
54959599516SKenneth E. Jansenc
55059599516SKenneth E. Jansen        if (Navier .eq. 1) then
55159599516SKenneth E. Jansenc
55259599516SKenneth E. Jansenc.... compute the inverse of deformation gradient
55359599516SKenneth E. Jansenc
55459599516SKenneth E. Jansen          dxidxb(:,1,1) =   dxdxib(:,2,2) * dxdxib(:,3,3)
55559599516SKenneth E. Jansen     &                    - dxdxib(:,3,2) * dxdxib(:,2,3)
55659599516SKenneth E. Jansen          dxidxb(:,1,2) =   dxdxib(:,3,2) * dxdxib(:,1,3)
55759599516SKenneth E. Jansen     &                    - dxdxib(:,1,2) * dxdxib(:,3,3)
55859599516SKenneth E. Jansen          dxidxb(:,1,3) =   dxdxib(:,1,2) * dxdxib(:,2,3)
55959599516SKenneth E. Jansen     &                    - dxdxib(:,1,3) * dxdxib(:,2,2)
56059599516SKenneth E. Jansen          temp          = one / ( dxidxb(:,1,1) * dxdxib(:,1,1)
56159599516SKenneth E. Jansen     &                          + dxidxb(:,1,2) * dxdxib(:,2,1)
56259599516SKenneth E. Jansen     &                          + dxidxb(:,1,3) * dxdxib(:,3,1) )
56359599516SKenneth E. Jansen          dxidxb(:,1,1) =  dxidxb(:,1,1) * temp
56459599516SKenneth E. Jansen          dxidxb(:,1,2) =  dxidxb(:,1,2) * temp
56559599516SKenneth E. Jansen          dxidxb(:,1,3) =  dxidxb(:,1,3) * temp
56659599516SKenneth E. Jansen          dxidxb(:,2,1) = (dxdxib(:,2,3) * dxdxib(:,3,1)
56759599516SKenneth E. Jansen     &                   - dxdxib(:,2,1) * dxdxib(:,3,3)) * temp
56859599516SKenneth E. Jansen          dxidxb(:,2,2) = (dxdxib(:,1,1) * dxdxib(:,3,3)
56959599516SKenneth E. Jansen     &                   - dxdxib(:,3,1) * dxdxib(:,1,3)) * temp
57059599516SKenneth E. Jansen          dxidxb(:,2,3) = (dxdxib(:,2,1) * dxdxib(:,1,3)
57159599516SKenneth E. Jansen     &                   - dxdxib(:,1,1) * dxdxib(:,2,3)) * temp
57259599516SKenneth E. Jansen          dxidxb(:,3,1) = (dxdxib(:,2,1) * dxdxib(:,3,2)
57359599516SKenneth E. Jansen     &                   - dxdxib(:,2,2) * dxdxib(:,3,1)) * temp
57459599516SKenneth E. Jansen          dxidxb(:,3,2) = (dxdxib(:,3,1) * dxdxib(:,1,2)
57559599516SKenneth E. Jansen     &                   - dxdxib(:,1,1) * dxdxib(:,3,2)) * temp
57659599516SKenneth E. Jansen          dxidxb(:,3,3) = (dxdxib(:,1,1) * dxdxib(:,2,2)
57759599516SKenneth E. Jansen     &                   - dxdxib(:,1,2) * dxdxib(:,2,1)) * temp
57859599516SKenneth E. Jansen
57959599516SKenneth E. Jansenc
58059599516SKenneth E. Jansenc.... compute local-grad-Y
58159599516SKenneth E. Jansenc
58259599516SKenneth E. Jansen          gl1yti = zero
58359599516SKenneth E. Jansen          gl2yti = zero
58459599516SKenneth E. Jansen          gl3yti = zero
58559599516SKenneth E. Jansenc
58659599516SKenneth E. Jansen          do n = 1, nshl
58759599516SKenneth E. Jansen            gl1yti(:) = gl1yti(:) + shglb(:,1,n) * ycl(:,n,id)
58859599516SKenneth E. Jansen            gl2yti(:) = gl2yti(:) + shglb(:,2,n) * ycl(:,n,id)
58959599516SKenneth E. Jansen            gl3yti(:) = gl3yti(:) + shglb(:,3,n) * ycl(:,n,id)
59059599516SKenneth E. Jansen          enddo
59159599516SKenneth E. Jansenc
59259599516SKenneth E. Jansenc.... convert local-grads to global-grads
59359599516SKenneth E. Jansenc
59459599516SKenneth E. Jansen          g1yti(:) = dxidxb(:,1,1) * gl1yti(:) +
59559599516SKenneth E. Jansen     &               dxidxb(:,2,1) * gl2yti(:) +
59659599516SKenneth E. Jansen     &               dxidxb(:,3,1) * gl3yti(:)
59759599516SKenneth E. Jansen          g2yti(:) = dxidxb(:,1,2) * gl1yti(:) +
59859599516SKenneth E. Jansen     &               dxidxb(:,2,2) * gl2yti(:) +
59959599516SKenneth E. Jansen     &               dxidxb(:,3,2) * gl3yti(:)
60059599516SKenneth E. Jansen          g3yti(:) = dxidxb(:,1,3) * gl1yti(:) +
60159599516SKenneth E. Jansen     &               dxidxb(:,2,3) * gl2yti(:) +
60259599516SKenneth E. Jansen     &               dxidxb(:,3,3) * gl3yti(:)
60359599516SKenneth E. Jansen
60459599516SKenneth E. Jansenc
60559599516SKenneth E. Jansenc.... end grad-Sclr
60659599516SKenneth E. Jansen        endif
60759599516SKenneth E. Jansenc
60859599516SKenneth E. Jansenc.... -------------------->  Boundary Conditions  <--------------------
60959599516SKenneth E. Jansenc
61059599516SKenneth E. Jansenc.... compute the Euler boundary conditions
61159599516SKenneth E. Jansenc
61259599516SKenneth E. Jansen        rou = zero
61359599516SKenneth E. Jansen        do n = 1, nshlb
61459599516SKenneth E. Jansen          nodlcl = lnode(n)
61559599516SKenneth E. Jansen          rou = rou + shpb(:,nodlcl) * BCB(:,n,1)
61659599516SKenneth E. Jansen        enddo
61759599516SKenneth E. Jansenc
61859599516SKenneth E. Jansenc.... impose scalar flux boundary conditions
61959599516SKenneth E. Jansen        SclrF = zero
62059599516SKenneth E. Jansen        do n=1,nshlb
62159599516SKenneth E. Jansen          nodlcl = lnode(n)
62259599516SKenneth E. Jansen          SclrF = SclrF + shpb(:,nodlcl) * BCB(:,n,ibb)
62359599516SKenneth E. Jansen        enddo
62459599516SKenneth E. Jansen
62559599516SKenneth E. Jansenc
62659599516SKenneth E. Jansenc.... flop count
62759599516SKenneth E. Jansenc
628*f4d0b58bSKenneth E. Jansen!      flops = flops + (27+18*nshl+14*nenbl)*npro
62959599516SKenneth E. Jansenc
63059599516SKenneth E. Jansenc.... return
63159599516SKenneth E. Jansenc
63259599516SKenneth E. Jansen        return
63359599516SKenneth E. Jansen        end
63459599516SKenneth E. Jansen
635