xref: /phasta/phSolver/compressible/rstatCheck.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine rstatCheck (res, ilwork,y,ac)
2*59599516SKenneth E. Jansenc
3*59599516SKenneth E. Jansenc----------------------------------------------------------------------
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansenc This subroutine calculates the statistics of the residual.
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc input:
8*59599516SKenneth E. Jansenc  res   (nshg,nflow)   : preconditioned residual
9*59599516SKenneth E. Jansenc
10*59599516SKenneth E. Jansenc output:
11*59599516SKenneth E. Jansenc  The time step, cpu-time and entropy-norm of the residual
12*59599516SKenneth E. Jansenc     are printed in the file HISTOR.DAT.
13*59599516SKenneth E. Jansenc
14*59599516SKenneth E. Jansenc
15*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
16*59599516SKenneth E. Jansenc----------------------------------------------------------------------
17*59599516SKenneth E. Jansenc
18*59599516SKenneth E. Jansen        include "common.h"
19*59599516SKenneth E. Jansen        include "mpif.h"
20*59599516SKenneth E. Jansen        include "auxmpi.h"
21*59599516SKenneth E. Jansenc
22*59599516SKenneth E. Jansen        dimension res(nshg,nflow)
23*59599516SKenneth E. Jansen        dimension rtmp(nshg), nrsmax(1), ilwork(nlwork)
24*59599516SKenneth E. Jansen        dimension Forin(4), Forout(4)
25*59599516SKenneth E. Jansen!SCATTER        dimension irecvcount(numpe), resvec(numpe)
26*59599516SKenneth E. Jansenc        integer TMRC
27*59599516SKenneth E. Jansen        real*8 y(nshg,ndof),ac(nshg,ndof)
28*59599516SKenneth E. Jansen        save ResLast
29*59599516SKenneth E. Jansen
30*59599516SKenneth E. Jansen        if (numpe == 1) nshgt=nshg   ! global = this processor
31*59599516SKenneth E. Jansenc
32*59599516SKenneth E. Jansenc.... ----------------------->  Convergence  <-------------------------
33*59599516SKenneth E. Jansenc
34*59599516SKenneth E. Jansenc.... compute the maximum residual and the corresponding node number
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansen        rtmp = zero
37*59599516SKenneth E. Jansen        do i = 1, nflow
38*59599516SKenneth E. Jansen          rtmp = rtmp + res(:,i)**2
39*59599516SKenneth E. Jansen        enddo
40*59599516SKenneth E. Jansen
41*59599516SKenneth E. Jansen        call sumgat (rtmp, 1, resnrm, ilwork)
42*59599516SKenneth E. Jansen
43*59599516SKenneth E. Jansen        resmaxl = maxval(rtmp)
44*59599516SKenneth E. Jansen
45*59599516SKenneth E. Jansen        irecvcount = 1
46*59599516SKenneth E. Jansen        resvec = resmaxl
47*59599516SKenneth E. Jansen        if (numpe > 1) then
48*59599516SKenneth E. Jansen        call MPI_ALLREDUCE (resvec, resmax, irecvcount,
49*59599516SKenneth E. Jansen     &                    MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD,
50*59599516SKenneth E. Jansen     &                    ierr)
51*59599516SKenneth E. Jansenc        call MPI_REDUCE_SCATTER (resvec, resmax, irecvcount,
52*59599516SKenneth E. Jansenc     &                    MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD,
53*59599516SKenneth E. Jansenc     &                    ierr)
54*59599516SKenneth E. Jansen        else
55*59599516SKenneth E. Jansen          resmax=resmaxl
56*59599516SKenneth E. Jansen        endif
57*59599516SKenneth E. Jansen        nrsmax = maxloc(rtmp)
58*59599516SKenneth E. Jansenc
59*59599516SKenneth E. Jansenc.... correct the residuals
60*59599516SKenneth E. Jansenc
61*59599516SKenneth E. Jansen        if (loctim(itseq) .eq. 0) then
62*59599516SKenneth E. Jansen          resnrm = resnrm
63*59599516SKenneth E. Jansen          resmax = resmax
64*59599516SKenneth E. Jansen        else
65*59599516SKenneth E. Jansen          resnrm = resnrm
66*59599516SKenneth E. Jansen          resmax = resmax
67*59599516SKenneth E. Jansen        endif
68*59599516SKenneth E. Jansenc
69*59599516SKenneth E. Jansenc.... approximate the number of entries
70*59599516SKenneth E. Jansenc
71*59599516SKenneth E. Jansen        totres = resnrm / float(nshgt)
72*59599516SKenneth E. Jansen        totres = sqrt(totres)
73*59599516SKenneth E. Jansen       if((iter.gt.1).and.(totres.gt.100.0*ResLast)) then !diverging
74*59599516SKenneth E. Jansen               call restar('out ',y,res) ! 'res' is used instead of 'ac'
75*59599516SKenneth E. Jansen               if(myrank.eq.0) write(*,*) 'ResLast totres', ResLast, totres
76*59599516SKenneth E. Jansen               if(myrank.eq.0) write(*,*) 'resmax', resmax
77*59599516SKenneth E. Jansen               call error('rstat    ','Diverge', iter)
78*59599516SKenneth E. Jansen       endif
79*59599516SKenneth E. Jansen       ResLast=totres
80*59599516SKenneth E. Jansen	ttim(68) = ttim(68) + secs(0.0)
81*59599516SKenneth E. Jansen
82*59599516SKenneth E. Jansenc
83*59599516SKenneth E. Jansenc.... return
84*59599516SKenneth E. Jansenc
85*59599516SKenneth E. Jansen        return
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansen        end
88*59599516SKenneth E. Jansen        subroutine rstatCheckSclr (rest, ilwork,y,ac)
89*59599516SKenneth E. Jansenc
90*59599516SKenneth E. Jansenc----------------------------------------------------------------------
91*59599516SKenneth E. Jansenc
92*59599516SKenneth E. Jansenc This subroutine calculates the statistics of the residual.
93*59599516SKenneth E. Jansenc
94*59599516SKenneth E. Jansenc input:
95*59599516SKenneth E. Jansenc  rest   (nshg)   : preconditioned residual
96*59599516SKenneth E. Jansenc
97*59599516SKenneth E. Jansenc output:
98*59599516SKenneth E. Jansenc  The time step, cpu-time and entropy-norm of the residual
99*59599516SKenneth E. Jansenc     are printed in the file HISTOR.DAT.
100*59599516SKenneth E. Jansenc
101*59599516SKenneth E. Jansenc
102*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
103*59599516SKenneth E. Jansenc----------------------------------------------------------------------
104*59599516SKenneth E. Jansenc
105*59599516SKenneth E. Jansen        include "common.h"
106*59599516SKenneth E. Jansen        include "mpif.h"
107*59599516SKenneth E. Jansen        include "auxmpi.h"
108*59599516SKenneth E. Jansenc
109*59599516SKenneth E. Jansen        dimension rest(nshg)
110*59599516SKenneth E. Jansen        dimension rtmp(nshg), nrsmax(1), ilwork(nlwork)
111*59599516SKenneth E. Jansen!SCATTER        dimension irecvcount(numpe), resvec(numpe)
112*59599516SKenneth E. Jansenc        integer TMRC
113*59599516SKenneth E. Jansen        real*8 y(nshg,ndof),ac(nshg,ndof)
114*59599516SKenneth E. Jansen        save ResLast
115*59599516SKenneth E. Jansen        save lstepLast
116*59599516SKenneth E. Jansen
117*59599516SKenneth E. Jansen	ttim(68) = ttim(68) - secs(0.0)
118*59599516SKenneth E. Jansen        if (numpe == 1) nshgt=nshg   ! global = this processor
119*59599516SKenneth E. Jansenc
120*59599516SKenneth E. Jansenc.... ----------------------->  Convergence  <-------------------------
121*59599516SKenneth E. Jansenc
122*59599516SKenneth E. Jansenc.... compute the maximum residual and the corresponding node number
123*59599516SKenneth E. Jansenc
124*59599516SKenneth E. Jansen        rtmp = zero
125*59599516SKenneth E. Jansen        rtmp = rtmp + rest**2
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansen        call sumgat (rtmp, 1, resnrm, ilwork)
128*59599516SKenneth E. Jansen
129*59599516SKenneth E. Jansen        resmaxl = maxval(rtmp)
130*59599516SKenneth E. Jansen
131*59599516SKenneth E. Jansencontinue on
132*59599516SKenneth E. Jansen
133*59599516SKenneth E. Jansen        irecvcount = 1
134*59599516SKenneth E. Jansen        resvec = resmaxl
135*59599516SKenneth E. Jansen        if (numpe > 1) then
136*59599516SKenneth E. Jansen        call MPI_ALLREDUCE (resvec, resmax, irecvcount,
137*59599516SKenneth E. Jansen     &                    MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD,
138*59599516SKenneth E. Jansen     &                    ierr)
139*59599516SKenneth E. Jansenc        call MPI_REDUCE_SCATTER (resvec, resmax, irecvcount,
140*59599516SKenneth E. Jansenc     &                    MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD,
141*59599516SKenneth E. Jansenc     &                    ierr)
142*59599516SKenneth E. Jansen        else
143*59599516SKenneth E. Jansen          resmax=resmaxl
144*59599516SKenneth E. Jansen        endif
145*59599516SKenneth E. Jansen        nrsmax = maxloc(rtmp)
146*59599516SKenneth E. Jansenc
147*59599516SKenneth E. Jansenc.... correct the residuals
148*59599516SKenneth E. Jansenc
149*59599516SKenneth E. Jansen        if (loctim(itseq) .eq. 0) then
150*59599516SKenneth E. Jansen          resnrm = resnrm
151*59599516SKenneth E. Jansen          resmax = resmax
152*59599516SKenneth E. Jansen        else
153*59599516SKenneth E. Jansen          resnrm = resnrm
154*59599516SKenneth E. Jansen          resmax = resmax
155*59599516SKenneth E. Jansen        endif
156*59599516SKenneth E. Jansenc
157*59599516SKenneth E. Jansenc.... approximate the number of entries
158*59599516SKenneth E. Jansenc
159*59599516SKenneth E. Jansen        totres = resnrm / float(nshgt)
160*59599516SKenneth E. Jansen        totres = sqrt(totres)
161*59599516SKenneth E. Jansen	if((lstep.gt.0).and.(lstepLast.eq.lstep)) then
162*59599516SKenneth E. Jansen           if(totres.gt.100.0*ResLast) then !diverging
163*59599516SKenneth E. Jansen               lstep = lstep+1
164*59599516SKenneth E. Jansen               ac(:,5) = rest(:) ! T dot in 'ac' is filled with scl. res
165*59599516SKenneth E. Jansen               call restar('out ',y,ac)
166*59599516SKenneth E. Jansen               if(myrank.eq.0) write(*,*) 'ResLast totres', ResLast, totres
167*59599516SKenneth E. Jansen               if(myrank.eq.0) write(*,*) 'resmax', resmax
168*59599516SKenneth E. Jansen               call error('rstatSclr','Diverge', iter)
169*59599516SKenneth E. Jansen           endif
170*59599516SKenneth E. Jansen	else
171*59599516SKenneth E. Jansen		lstepLast=lstep
172*59599516SKenneth E. Jansen	endif
173*59599516SKenneth E. Jansen       ResLast=totres
174*59599516SKenneth E. Jansenc
175*59599516SKenneth E. Jansenc.... return
176*59599516SKenneth E. Jansenc
177*59599516SKenneth E. Jansen        return
178*59599516SKenneth E. Jansenc
179*59599516SKenneth E. Jansen        end
180