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