xref: /phasta/phSolver/compressible/asigmr.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08) !
1        subroutine AsIGMR (y,       ac,      x,        xmudmi,
2     &                     shp,     shgl,    ien,
3     &                     mater,   res,     rmes,
4     &                     BDiag,   qres,    EGmass,   rerr)
5c
6c----------------------------------------------------------------------
7c
8c This routine computes and assembles the data corresponding to the
9c  interior elements.
10c
11c Zdenek Johan, Winter 1991.  (Fortran 90)
12c----------------------------------------------------------------------
13c
14        use rlssave     ! Use the resolved Leonard stresses at the nodes.
15        use timedataC    ! time series
16        use specialBC    ! get ytarget to localize and send down
17        include "common.h"
18c
19        dimension y(nshg,ndofl),            ac(nshg,ndofl),
20     &            x(numnp,nsd),
21     &            shp(nshl,MAXQPT),
22     &            shgl(nsd,nshl,MAXQPT),
23     &            ien(npro,nshl),
24     &            mater(npro),              res(nshg,nflow),
25     &            rmes(nshg,nflow),         BDiag(nshg,nflow,nflow),
26     &            qres(nshg,idflx)
27
28c
29        dimension ycl(npro,nshl,ndofl),     acl(npro,nshl,ndof),
30     &            xl(npro,nenl,nsd),        ytargetl(npro,nshl,nflow),
31     &            rl(npro,nshl,nflow),      rml(npro,nshl,nflow),
32     &            BDiagl(npro,nshl,nflow,nflow),
33     &            ql(npro,nshl,idflx)
34c
35        dimension  xmudmi(npro,ngauss)
36        dimension sgn(npro,nshl),  EGmass(npro,nedof,nedof)
37
38        dimension rlsl(npro,nshl,6)
39        real*8 rerrl(npro,nshl,6), rerr(nshg,10)
40
41c
42c.... create the matrix of mode signs for the hierarchic basis
43c     functions.
44c
45c
46        if (ipord .gt. 1) then
47           call getsgn(ien,sgn)
48        endif
49c
50c.... gather the variables
51c
52        call localy(y,      ycl,     ien,    ndofl,  'gather  ')
53        call localy(ac,    acl,     ien,    ndofl,  'gather  ')
54        call localx(x,      xl,     ien,    nsd,    'gather  ')
55        call local (qres,   ql,     ien,    idflx,  'gather  ')
56
57        if(matflg(5,1).ge.4 )
58     &   call localy (ytarget,   ytargetl,  ien,   nflow,  'gather  ')
59
60
61        if( (iLES.gt.10).and.(iLES.lt.20)) then  ! bardina
62           call local (rls, rlsl,     ien,       6, 'gather  ')
63        else
64           rlsl = zero
65        endif
66c
67c.... get the element residuals, LHS matrix, and preconditioner
68c
69        rl     = zero
70        BDiagl = zero
71
72        if(ierrcalc.eq.1) rerrl = zero
73        ttim(31) = ttim(31) - secs(0.0)
74
75        call e3  (ycl,     ycl,     acl,     shp,
76     &            shgl,    xl,      rl,      rml,   xmudmi,
77     &            BDiagl,  ql,      sgn,     rlsl,  EGmass,
78     &            rerrl,   ytargetl)
79
80        ttim(31) = ttim(31) + secs(0.0)
81c
82c.... assemble the residual and modified residual
83c
84        call local (res,    rl,     ien,    nflow,  'scatter ')
85c
86        if ( ierrcalc .eq. 1 ) then
87           call local (rerr, rerrl,  ien, 6, 'scatter ')
88        endif
89c
90c.... extract and assemble the Block-Diagonal (see note in elmgmr, line 280)
91c
92        if (iprec .ne. 0) then
93           do i = 1, nshl
94              do j = 1, nflow
95                 i0 = (i - 1) * nflow + j
96                 do k = 1, nflow
97                    j0 = (i - 1) * nflow + k
98                    BDiagl(:,i,j,k) = EGmass(:,i0,j0)
99                 enddo
100              enddo
101           enddo
102           call local (BDiag,  BDiagl, ien, nflow*nflow, 'scatter ')
103        endif
104
105c
106c... call timeseries
107c
108
109        if (exts) then
110           if ((iter.eq.1).and.(mod(lstep,freq).eq.0)) then
111              call timeseries(ycl,xl,ien,sgn)
112           endif
113        endif
114
115c
116c.... end
117c
118        return
119        end
120c
121c
122c
123        subroutine AsIGMRSclr (y,       ac,
124     &                         x,       elDwl,
125     &                         shp,     shgl,      ien,
126     &                         mater,   rest,      rmest,
127     &                         qrest,   EGmasst,   Diag)
128c
129c----------------------------------------------------------------------
130c
131c This routine computes and assembles the data corresponding to the
132c  interior elements.
133c
134c Zdenek Johan, Winter 1991.  (Fortran 90)
135c----------------------------------------------------------------------
136c
137        use turbSA
138        include "common.h"
139c
140        dimension y(nshg,ndof),
141     &            ac(nshg,ndof),
142     &            x(numnp,nsd),
143     &            shp(nshl,MAXQPT),        shgl(nsd,nshl,MAXQPT),
144     &            ien(npro,nshl),
145     &            mater(npro),            rest(nshg),
146     &            rmest(nshg),            Diag(nshg),
147     &            qrest(nshg)
148
149c
150        dimension ycl(npro,nshl,ndof),
151     &            acl(npro,nshl,ndof),    dwl(npro,nenl),
152     &            xl(npro,nenl,nsd),      Diagl(npro,nshl),
153     &            rtl(npro,nshl),         rmtl(npro,nshl),
154     &            qtl(npro,nshl),         sgn(npro,nshl)
155c
156        dimension EGmasst(npro,nshape, nshape)
157        real*8    elDwl(npro)
158c.... create the matrix of mode signs for the hierarchic basis
159c     functions.
160c
161        call getsgn(ien,sgn)
162c
163c
164c.... gather the variables
165c
166        call localy (y,       ycl,      ien,    ndof,  'gather  ')
167        call localy (ac,      acl,     ien,    ndof,  'gather  ')
168        call localx (x,       xl,      ien,    nsd,   'gather  ')
169c       call local (qrest,   qtl,     ien,    1,     'gather  ')
170        if (iRANS .lt. 0) then
171           call localx (d2wall,   dwl,     ien,    1,     'gather  ')
172        endif
173c
174c.... get the element residuals, LHS matrix, and preconditioner
175c
176        rtl     = zero
177        Diagl   = zero
178
179        ttim(31) = ttim(31) - tmr()
180
181        call e3Sclr (ycl,     acl,
182     &               dwl,    elDwl,   shp,
183     &               sgn,    shgl,    xl,
184     &               rtl,    rmtl,
185     &               qtl,    EGmasst )
186
187        ttim(31) = ttim(31) + tmr()
188c
189c.... assemble the residual and modified residual
190c
191        call local (rest,    rtl,     ien,    1,  'scatter ')
192c
193c.... extract and assemble the Diagonal
194c
195        if (iprec .ne. 0) then
196           do i=1,nshl
197              Diagl(:,i)=EGmassT(:,i,i)
198           enddo
199           call local(Diag, Diagl, ien, 1, 'scatter ')
200        endif
201c
202c.... end
203c
204        return
205        end
206
207
208