xref: /phasta/phSolver/incompressible/asiq.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine AsIq (y,       x,       shp,
2*59599516SKenneth E. Jansen     &                   shgl,    ien,     xmudmi,
3*59599516SKenneth E. Jansen     &                   qres,    rmass    )
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansenc----------------------------------------------------------------------
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc This routine computes and assembles the data corresponding to the
8*59599516SKenneth E. Jansenc interior elements for the global reconstruction of the diffusive
9*59599516SKenneth E. Jansenc flux vector.
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansenc input:
12*59599516SKenneth E. Jansenc     y     (numnp,ndof)        : Y variables
13*59599516SKenneth E. Jansenc     x     (numnp,nsd)         : nodal coordinates
14*59599516SKenneth E. Jansenc     shp   (nen,nintg)         : element shape-functions
15*59599516SKenneth E. Jansenc     shgl  (nsd,nen,nintg)     : element local shape-function gradients
16*59599516SKenneth E. Jansenc     ien   (npro)              : nodal connectivity array
17*59599516SKenneth E. Jansenc
18*59599516SKenneth E. Jansenc output:
19*59599516SKenneth E. Jansenc     qres  (numnp,nsd,nsd)  : residual vector for diffusive flux
20*59599516SKenneth E. Jansenc     rmass  (numnp)            : lumped mass matrix
21*59599516SKenneth E. Jansenc
22*59599516SKenneth E. Jansenc----------------------------------------------------------------------
23*59599516SKenneth E. Jansenc
24*59599516SKenneth E. Jansen        use turbsa      ! access to d2wall
25*59599516SKenneth E. Jansen        include "common.h"
26*59599516SKenneth E. Jansenc
27*59599516SKenneth E. Jansen        dimension y(nshg,ndof),               x(numnp,nsd),
28*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),         shgl(nsd,nshl,ngauss),
29*59599516SKenneth E. Jansen     &            ien(npro,nshl),      dwl(npro,nenl),
30*59599516SKenneth E. Jansen     &            qres(nshg,idflx),    rmass(nshg)
31*59599516SKenneth E. Jansenc
32*59599516SKenneth E. Jansen        dimension yl(npro,nshl,ndof),          xl(npro,nenl,nsd),
33*59599516SKenneth E. Jansen     &            ql(npro,nshl,idflx),  rmassl(npro,nshl),
34*59599516SKenneth E. Jansen     &            xmudmi(npro,ngauss)
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansen        dimension sgn(npro,nshl)
37*59599516SKenneth E. Jansenc
38*59599516SKenneth E. Jansenc.... create the matrix of mode signs for the hierarchic basis
39*59599516SKenneth E. Jansenc     functions.
40*59599516SKenneth E. Jansenc
41*59599516SKenneth E. Jansen        do i=1,nshl
42*59599516SKenneth E. Jansen           where ( ien(:,i) < 0 )
43*59599516SKenneth E. Jansen              sgn(:,i) = -one
44*59599516SKenneth E. Jansen           elsewhere
45*59599516SKenneth E. Jansen              sgn(:,i) = one
46*59599516SKenneth E. Jansen           endwhere
47*59599516SKenneth E. Jansen        enddo
48*59599516SKenneth E. Jansen
49*59599516SKenneth E. Jansenc
50*59599516SKenneth E. Jansenc.... gather the variables
51*59599516SKenneth E. Jansenc
52*59599516SKenneth E. Jansen
53*59599516SKenneth E. Jansen        call localy(y,      yl,     ien,    ndof,   'gather  ')
54*59599516SKenneth E. Jansen        call localx (x,      xl,     ien,    nsd,    'gather  ')
55*59599516SKenneth E. Jansen        if (iRANS .eq. -2) then ! kay-epsilon
56*59599516SKenneth E. Jansen           call localx (d2wall,   dwl,     ien,    1,     'gather  ')
57*59599516SKenneth E. Jansen        endif
58*59599516SKenneth E. Jansenc
59*59599516SKenneth E. Jansenc.... get the element residuals
60*59599516SKenneth E. Jansenc
61*59599516SKenneth E. Jansen        ql     = zero
62*59599516SKenneth E. Jansen        rmassl = zero
63*59599516SKenneth E. Jansen
64*59599516SKenneth E. Jansen        call e3q  (yl,         dwl,      shp,      shgl,
65*59599516SKenneth E. Jansen     &             xl,         ql,       rmassl,
66*59599516SKenneth E. Jansen     &             xmudmi,     sgn  )
67*59599516SKenneth E. Jansen
68*59599516SKenneth E. Jansenc
69*59599516SKenneth E. Jansenc.... assemble the diffusive flux residual
70*59599516SKenneth E. Jansenc
71*59599516SKenneth E. Jansen        call local (qres,   ql,  ien,  idflx,  'scatter ')
72*59599516SKenneth E. Jansen        call local (rmass,  rmassl, ien,  1,          'scatter ')
73*59599516SKenneth E. Jansenc
74*59599516SKenneth E. Jansenc.... end
75*59599516SKenneth E. Jansenc
76*59599516SKenneth E. Jansen        return
77*59599516SKenneth E. Jansen        end
78*59599516SKenneth E. Jansen
79*59599516SKenneth E. Jansen
80*59599516SKenneth E. Jansenc
81*59599516SKenneth E. Jansenc----------------------------------------------------------------------
82*59599516SKenneth E. Jansenc
83*59599516SKenneth E. Jansenc This routine computes and assembles the data corresponding to the
84*59599516SKenneth E. Jansenc interior elements for the global reconstruction of the diffusive
85*59599516SKenneth E. Jansenc flux vector.
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansenc----------------------------------------------------------------------
88*59599516SKenneth E. Jansen        subroutine AsIqSclr (y,       x,       shp,
89*59599516SKenneth E. Jansen     &                       shgl,    ien,     qres,
90*59599516SKenneth E. Jansen     &                       rmass    )
91*59599516SKenneth E. Jansenc
92*59599516SKenneth E. Jansen        use turbsa      ! access to d2wall
93*59599516SKenneth E. Jansen        include "common.h"
94*59599516SKenneth E. Jansenc
95*59599516SKenneth E. Jansen        dimension y(nshg,ndof),             x(numnp,nsd),
96*59599516SKenneth E. Jansen     &            shp(nshl,ngauss),         shgl(nsd,nshl,ngauss),
97*59599516SKenneth E. Jansen     &            ien(npro,nshl),      dwl(npro,nenl),
98*59599516SKenneth E. Jansen     &            qres(nshg,nsd),           rmass(nshg)
99*59599516SKenneth E. Jansenc
100*59599516SKenneth E. Jansen        dimension yl(npro,nshl,ndof),       xl(npro,nenl,nsd),
101*59599516SKenneth E. Jansen     &            ql(npro,nshl,nsd),        rmassl(npro,nshl)
102*59599516SKenneth E. Jansenc
103*59599516SKenneth E. Jansen        dimension sgn(npro,nshl)
104*59599516SKenneth E. Jansen
105*59599516SKenneth E. Jansen        if (ipord .gt. 1) then
106*59599516SKenneth E. Jansen           call getsgn(ien,sgn)
107*59599516SKenneth E. Jansen        endif
108*59599516SKenneth E. Jansenc
109*59599516SKenneth E. Jansenc.... gather the variables
110*59599516SKenneth E. Jansenc
111*59599516SKenneth E. Jansen        call localy(y,      yl,     ien,    ndof,   'gather  ')
112*59599516SKenneth E. Jansen        call localx (x,      xl,     ien,    nsd,    'gather  ')
113*59599516SKenneth E. Jansen        if (iRANS .eq. -2) then ! kay-epsilon
114*59599516SKenneth E. Jansen           call localx (d2wall,   dwl,     ien,    1,     'gather  ')
115*59599516SKenneth E. Jansen        endif
116*59599516SKenneth E. Jansenc
117*59599516SKenneth E. Jansenc.... get the element residuals
118*59599516SKenneth E. Jansenc
119*59599516SKenneth E. Jansen        ql     = zero
120*59599516SKenneth E. Jansen        rmassl = zero
121*59599516SKenneth E. Jansen
122*59599516SKenneth E. Jansen        call e3qSclr  (yl,      dwl,    shp,    shgl,
123*59599516SKenneth E. Jansen     &                 xl,      ql,     rmassl,
124*59599516SKenneth E. Jansen     &                 sgn             )
125*59599516SKenneth E. Jansen
126*59599516SKenneth E. Jansenc
127*59599516SKenneth E. Jansenc.... assemble the temperature diffusive flux residual
128*59599516SKenneth E. Jansenc
129*59599516SKenneth E. Jansen        call local (qres,   ql,  ien,  nsd,  'scatter ')
130*59599516SKenneth E. Jansen        call local (rmass,  rmassl, ien,  1, 'scatter ')
131*59599516SKenneth E. Jansenc
132*59599516SKenneth E. Jansenc.... end
133*59599516SKenneth E. Jansenc
134*59599516SKenneth E. Jansen        return
135*59599516SKenneth E. Jansen        end
136*59599516SKenneth E. Jansen
137