xref: /phasta/phSolver/common/errorsmooth.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
1*10167291SKenneth E. Jansen      subroutine errsmooth(rerr,   x,     iper,   ilwork,
2*10167291SKenneth E. Jansen     &                     shp,    shgl,  iBC)
3*10167291SKenneth E. Jansenc
4*10167291SKenneth E. Jansen        use pointer_data
5*10167291SKenneth E. Jansenc
6*10167291SKenneth E. Jansen        include "common.h"
7*10167291SKenneth E. Jansen        include "mpif.h"
8*10167291SKenneth E. Jansenc
9*10167291SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
10*10167291SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
11*10167291SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
12*10167291SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
13*10167291SKenneth E. Jansenc
14*10167291SKenneth E. Jansen        dimension rerrsm(nshg, 10), rerr(nshg,10), rmass(nshg),
15*10167291SKenneth E. Jansen     &            x(nshg,3)
16*10167291SKenneth E. Jansenc
17*10167291SKenneth E. Jansen        dimension ilwork(nlwork), iBC(nshg), iper(nshg)
18*10167291SKenneth E. Jansen
19*10167291SKenneth E. Jansen        real*8, allocatable :: tmpshp(:,:), tmpshgl(:,:,:)
20*10167291SKenneth E. Jansen        real*8, allocatable :: tmpshpb(:,:), tmpshglb(:,:,:)
21*10167291SKenneth E. Jansen
22*10167291SKenneth E. Jansenc
23*10167291SKenneth E. Jansenc loop over element blocks for the global reconstruction
24*10167291SKenneth E. Jansenc of the smoothed error and lumped mass matrix, rmass
25*10167291SKenneth E. Jansenc
26*10167291SKenneth E. Jansen        rerrsm = zero
27*10167291SKenneth E. Jansen        rmass = zero
28*10167291SKenneth E. Jansen
29*10167291SKenneth E. Jansen        do iblk = 1, nelblk
30*10167291SKenneth E. Jansenc
31*10167291SKenneth E. Jansenc.... set up the parameters
32*10167291SKenneth E. Jansenc
33*10167291SKenneth E. Jansen          nenl   = lcblk(5,iblk)   ! no. of vertices per element
34*10167291SKenneth E. Jansen          iel    = lcblk(1,iblk)
35*10167291SKenneth E. Jansen          lelCat = lcblk(2,iblk)
36*10167291SKenneth E. Jansen          lcsyst = lcblk(3,iblk)
37*10167291SKenneth E. Jansen          iorder = lcblk(4,iblk)
38*10167291SKenneth E. Jansen          nenl   = lcblk(5,iblk)   ! no. of vertices per element
39*10167291SKenneth E. Jansen          nshl   = lcblk(10,iblk)
40*10167291SKenneth E. Jansen          mattyp = lcblk(7,iblk)
41*10167291SKenneth E. Jansen          ndofl  = lcblk(8,iblk)
42*10167291SKenneth E. Jansen          nsymdl = lcblk(9,iblk)
43*10167291SKenneth E. Jansen          npro   = lcblk(1,iblk+1) - iel
44*10167291SKenneth E. Jansen          ngauss = nint(lcsyst)
45*10167291SKenneth E. Jansenc
46*10167291SKenneth E. Jansenc.... compute and assemble diffusive flux vector residual, qres,
47*10167291SKenneth E. Jansenc     and lumped mass matrix, rmass
48*10167291SKenneth E. Jansen
49*10167291SKenneth E. Jansen          allocate (tmpshp(nshl,MAXQPT))
50*10167291SKenneth E. Jansen          allocate (tmpshgl(nsd,nshl,MAXQPT))
51*10167291SKenneth E. Jansen
52*10167291SKenneth E. Jansen          tmpshp(1:nshl,:) = shp(lcsyst,1:nshl,:)
53*10167291SKenneth E. Jansen          tmpshgl(:,1:nshl,:) = shgl(lcsyst,:,1:nshl,:)
54*10167291SKenneth E. Jansen
55*10167291SKenneth E. Jansen          call smooth (rerr,                x,
56*10167291SKenneth E. Jansen     &               tmpshp,
57*10167291SKenneth E. Jansen     &               tmpshgl,
58*10167291SKenneth E. Jansen     &               mien(iblk)%p,
59*10167291SKenneth E. Jansen     &               rerrsm,
60*10167291SKenneth E. Jansen     &               rmass)
61*10167291SKenneth E. Jansen
62*10167291SKenneth E. Jansen          deallocate ( tmpshp )
63*10167291SKenneth E. Jansen          deallocate ( tmpshgl )
64*10167291SKenneth E. Jansen       enddo
65*10167291SKenneth E. Jansenc
66*10167291SKenneth E. Jansen       if (numpe > 1) then
67*10167291SKenneth E. Jansen          call commu (rerrsm , ilwork,  10   , 'in ')
68*10167291SKenneth E. Jansen          call commu (rmass  , ilwork,  1    , 'in ')
69*10167291SKenneth E. Jansen       endif
70*10167291SKenneth E. Jansenc
71*10167291SKenneth E. Jansenc.... take care of periodic boundary conditions
72*10167291SKenneth E. Jansenc
73*10167291SKenneth E. Jansen        do j= 1,nshg
74*10167291SKenneth E. Jansen          if ((btest(iBC(j),10))) then
75*10167291SKenneth E. Jansen            i = iper(j)
76*10167291SKenneth E. Jansen            rmass(i) = rmass(i) + rmass(j)
77*10167291SKenneth E. Jansen            rerrsm(i,:) = rerrsm(i,:) + rerrsm(j,:)
78*10167291SKenneth E. Jansen          endif
79*10167291SKenneth E. Jansen        enddo
80*10167291SKenneth E. Jansen
81*10167291SKenneth E. Jansen        do j= 1,nshg
82*10167291SKenneth E. Jansen          if ((btest(iBC(j),10))) then
83*10167291SKenneth E. Jansen            i = iper(j)
84*10167291SKenneth E. Jansen            rmass(j) = rmass(i)
85*10167291SKenneth E. Jansen            rerrsm(j,:) = rerrsm(i,:)
86*10167291SKenneth E. Jansen          endif
87*10167291SKenneth E. Jansen        enddo
88*10167291SKenneth E. Jansenc
89*10167291SKenneth E. Jansenc.... invert the diagonal mass matrix and find q
90*10167291SKenneth E. Jansenc
91*10167291SKenneth E. Jansen        rmass = one/rmass
92*10167291SKenneth E. Jansen
93*10167291SKenneth E. Jansen       do i=1, 10
94*10167291SKenneth E. Jansen          rerrsm(:,i) = rmass*rerrsm(:,i)
95*10167291SKenneth E. Jansen       enddo
96*10167291SKenneth E. Jansen       if(numpe > 1) then
97*10167291SKenneth E. Jansen          call commu (rerrsm, ilwork, 10, 'out')
98*10167291SKenneth E. Jansen       endif
99*10167291SKenneth E. Jansenc
100*10167291SKenneth E. Jansenc      copy the smoothed error overwriting the original error.
101*10167291SKenneth E. Jansenc
102*10167291SKenneth E. Jansen
103*10167291SKenneth E. Jansen       rerr = rerrsm
104*10167291SKenneth E. Jansen
105*10167291SKenneth E. Jansen       return
106*10167291SKenneth E. Jansen       end
107*10167291SKenneth E. Jansen
108*10167291SKenneth E. Jansen        subroutine smooth (rerr,       x,       shp,
109*10167291SKenneth E. Jansen     &                     shgl,       ien,
110*10167291SKenneth E. Jansen     &                     rerrsm,     rmass    )
111*10167291SKenneth E. Jansenc
112*10167291SKenneth E. Jansenc----------------------------------------------------------------------
113*10167291SKenneth E. Jansenc
114*10167291SKenneth E. Jansenc This routine computes and assembles the data corresponding to the
115*10167291SKenneth E. Jansenc interior elements for the global reconstruction of the diffusive
116*10167291SKenneth E. Jansenc flux vector.
117*10167291SKenneth E. Jansenc
118*10167291SKenneth E. Jansenc input:
119*10167291SKenneth E. Jansenc     y     (nshg,ndof)        : Y variables
120*10167291SKenneth E. Jansenc     x     (numnp,nsd)         : nodal coordinates
121*10167291SKenneth E. Jansenc     shp   (nshape,ngauss)     : element shape-functions
122*10167291SKenneth E. Jansenc     shgl  (nsd,nshape,ngauss) : element local shape-function gradients
123*10167291SKenneth E. Jansenc     ien   (npro)              : nodal connectivity array
124*10167291SKenneth E. Jansenc
125*10167291SKenneth E. Jansenc output:
126*10167291SKenneth E. Jansenc     qres  (nshg,nflow-1,nsd)  : residual vector for diffusive flux
127*10167291SKenneth E. Jansenc     rmass  (nshg)            : lumped mass matrix
128*10167291SKenneth E. Jansenc
129*10167291SKenneth E. Jansenc----------------------------------------------------------------------
130*10167291SKenneth E. Jansenc
131*10167291SKenneth E. Jansen        include "common.h"
132*10167291SKenneth E. Jansenc
133*10167291SKenneth E. Jansen        dimension rerr(nshg,10),               x(numnp,nsd),
134*10167291SKenneth E. Jansen     &            shp(nshl,maxsh),
135*10167291SKenneth E. Jansen     &            shgl(nsd,nshl,maxsh),
136*10167291SKenneth E. Jansen     &            ien(npro,nshl),
137*10167291SKenneth E. Jansen     &            rerrsm(nshg,10),    rmass(nshg)
138*10167291SKenneth E. Jansenc
139*10167291SKenneth E. Jansenc.... element level declarations
140*10167291SKenneth E. Jansenc
141*10167291SKenneth E. Jansen        dimension rerrl(npro,nshl,10),        xl(npro,nenl,nsd),
142*10167291SKenneth E. Jansen     &            rerrsml(npro,nshl,10),       rmassl(npro,nshl)
143*10167291SKenneth E. Jansenc
144*10167291SKenneth E. Jansen        dimension sgn(npro,nshl),          shape(npro,nshl),
145*10167291SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),    WdetJ(npro),
146*10167291SKenneth E. Jansen     &            dxidx(npro,nsd,nsd),     shg(npro,nshl,nsd)
147*10167291SKenneth E. Jansenc
148*10167291SKenneth E. Jansen        dimension error(npro,10)
149*10167291SKenneth E. Jansenc
150*10167291SKenneth E. Jansenc.... create the matrix of mode signs for the hierarchic basis
151*10167291SKenneth E. Jansenc     functions.
152*10167291SKenneth E. Jansenc
153*10167291SKenneth E. Jansen        if (ipord .gt. 1) then
154*10167291SKenneth E. Jansen           call getsgn(ien,sgn)
155*10167291SKenneth E. Jansen        endif
156*10167291SKenneth E. Jansenc
157*10167291SKenneth E. Jansenc.... gather the variables
158*10167291SKenneth E. Jansenc
159*10167291SKenneth E. Jansen
160*10167291SKenneth E. Jansen        call local(rerr,   rerrl,  ien,    10,   'gather  ')
161*10167291SKenneth E. Jansen        call localx(x,      xl,     ien,    nsd,    'gather  ')
162*10167291SKenneth E. Jansenc
163*10167291SKenneth E. Jansenc.... get the element residuals
164*10167291SKenneth E. Jansenc
165*10167291SKenneth E. Jansen        rerrsml     = zero
166*10167291SKenneth E. Jansen        rmassl      = zero
167*10167291SKenneth E. Jansen
168*10167291SKenneth E. Jansenc
169*10167291SKenneth E. Jansenc.... loop through the integration points
170*10167291SKenneth E. Jansenc
171*10167291SKenneth E. Jansen
172*10167291SKenneth E. Jansen
173*10167291SKenneth E. Jansen        do intp = 1, ngauss
174*10167291SKenneth E. Jansen        if (Qwt(lcsyst,intp) .eq. zero) cycle          ! precaution
175*10167291SKenneth E. Jansenc
176*10167291SKenneth E. Jansenc.... create a matrix of shape functions (and derivatives) for each
177*10167291SKenneth E. Jansenc     element at this quadrature point. These arrays will contain
178*10167291SKenneth E. Jansenc     the correct signs for the hierarchic basis
179*10167291SKenneth E. Jansenc
180*10167291SKenneth E. Jansen        call getshp(shp,          shgl,      sgn,
181*10167291SKenneth E. Jansen     &              shape,        shdrv)
182*10167291SKenneth E. Jansenc
183*10167291SKenneth E. Jansen        call e3metric( xl,         shdrv,        dxidx,
184*10167291SKenneth E. Jansen     &                 shg,        WdetJ)
185*10167291SKenneth E. Jansen        error=zero
186*10167291SKenneth E. Jansen        do n = 1, nshl
187*10167291SKenneth E. Jansen           do i=1,10
188*10167291SKenneth E. Jansen              error(:,i)=error(:,i) + shape(:,n) * rerrl(:,n,i)
189*10167291SKenneth E. Jansen           enddo
190*10167291SKenneth E. Jansen        enddo
191*10167291SKenneth E. Jansen        do i=1,nshl
192*10167291SKenneth E. Jansen           do j=1,10
193*10167291SKenneth E. Jansen              rerrsml(:,i,j)  = rerrsml(:,i,j)
194*10167291SKenneth E. Jansen     &                       + shape(:,i)*WdetJ*error(:,j)
195*10167291SKenneth E. Jansen           enddo
196*10167291SKenneth E. Jansen
197*10167291SKenneth E. Jansen           rmassl(:,i) = rmassl(:,i) + shape(:,i)*WdetJ
198*10167291SKenneth E. Jansen        enddo
199*10167291SKenneth E. Jansen
200*10167291SKenneth E. Jansenc.... end of the loop over integration points
201*10167291SKenneth E. Jansenc
202*10167291SKenneth E. Jansen      enddo
203*10167291SKenneth E. Jansenc
204*10167291SKenneth E. Jansenc.... assemble the diffusive flux residual
205*10167291SKenneth E. Jansenc
206*10167291SKenneth E. Jansen        call local (rerrsm,   rerrsml,  ien,  10,'scatter ')
207*10167291SKenneth E. Jansen        call local (rmass,   rmassl,  ien,  1,  'scatter ')
208*10167291SKenneth E. Jansenc
209*10167291SKenneth E. Jansen
210*10167291SKenneth E. Jansen      return
211*10167291SKenneth E. Jansen      end
212