/[MITgcm]/MITgcm/pkg/grdchk/grdchk_print.F
ViewVC logotype

Diff of /MITgcm/pkg/grdchk/grdchk_print.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.1 by heimbach, Fri Jul 13 13:08:17 2001 UTC revision 1.2 by heimbach, Fri Jul 13 14:50:46 2001 UTC
# Line 0  Line 1 
1    C $Header$
2    
3    #include "CPP_OPTIONS.h"
4    
5    
6          subroutine grdchk_print(
7         I                         ichknum,
8         I                         ierr_grdchk,
9         I                         mythid
10         &                       )
11    
12    c     ==================================================================
13    c     SUBROUTINE grdchk_print
14    c     ==================================================================
15    c
16    c     o Print the results of the gradient check.
17    c
18    c     started: Christian Eckert eckert@mit.edu 08-Mar-2000
19    c     continued: heimbach@mit.edu: 13-Jun-2001
20    c
21    c     ==================================================================
22    c     SUBROUTINE grdchk_print
23    c     ==================================================================
24    
25          implicit none
26    
27    c     == global variables ==
28    
29    #include "EEPARAMS.h"
30    #include "SIZE.h"
31    #include "GRID.h"
32    #include "grdchk.h"
33    
34    c     == routine arguments ==
35    
36          integer ichknum
37          integer ierr_grdchk
38          integer mythid
39    
40    #ifdef ALLOW_GRADIENT_CHECK
41    c     == local variables ==
42    
43          _RL fcref
44          _RL fcpert
45          _RL xxmemo_ref
46          _RL xxmemo_pert
47          _RL gfd
48          _RL adxxmemo
49          _RL ratio
50    
51          integer i
52          integer itile
53          integer jtile
54          integer itilepos
55          integer jtilepos
56          integer layer
57          integer icomp
58          integer ierr
59    
60          integer numchecks
61    
62          character*(max_len_mbuf) msgbuf
63    
64    c     == end of interface ==
65    
66    c--   Print header.
67          write(msgbuf,'(a)')
68         &' '
69          call print_message( msgbuf, standardmessageunit,
70         &                    SQUEEZE_RIGHT , mythid)
71          write(msgbuf,'(a)')
72         &'// ======================================================='
73          call print_message( msgbuf, standardmessageunit,
74         &                    SQUEEZE_RIGHT , mythid)
75          write(msgbuf,'(a)')
76         &'// Gradient check results  >>> START <<<'
77          call print_message( msgbuf, standardmessageunit,
78         &                    SQUEEZE_RIGHT , mythid)
79          write(msgbuf,'(a)')
80         &'// ======================================================='
81          call print_message( msgbuf, standardmessageunit,
82         &                    SQUEEZE_RIGHT , mythid)
83          write(msgbuf,'(a)')
84         &' '
85          call print_message( msgbuf, standardmessageunit,
86         &                    SQUEEZE_RIGHT , mythid)
87          write(msgbuf,'(a)')
88         &' '
89          call print_message( msgbuf, standardmessageunit,
90         &                    SQUEEZE_RIGHT , mythid)
91    
92          write(msgbuf,'(a,e10.3)')
93         &' EPS = ',grdchk_eps
94          call print_message( msgbuf, standardmessageunit,
95         &                    SQUEEZE_RIGHT , mythid)
96    
97          write(msgbuf,'(a,a)')
98         &'     I         X(I)           FC          FC1-FC2/EPS ',
99         &'     GRAD(FC)   DIFFERENCE/RATIO '
100          call print_message( msgbuf, standardmessageunit,
101         &                    SQUEEZE_RIGHT , mythid)
102    
103    c--   Individual checks.
104          if ( ierr_grdchk .eq. 0 ) then
105             numchecks = ichknum - 1
106          else
107             numchecks = maxgrdchecks
108          endif
109    
110          do i = 1, numchecks
111            xxmemo_ref   = xxmemref  (i)
112            xxmemo_pert  = xxmempert (i)
113            adxxmemo     = adxxmem   (i)
114            fcref    = fcrmem  (i)
115            fcpert   = fcpmem  (i)
116            gfd      = gfdmem  (i)
117            ratio    = ratiomem(i)
118            itile    = bimem   (i)
119            jtile    = bjmem   (i)
120            itilepos = ilocmem (i)
121            jtilepos = jlocmem (i)
122            layer    = klocmem (i)
123            icomp    = icompmem(i)
124            ierr     = ierrmem (i)
125    
126            write(msgbuf,'(A,4(I8),2(1x,D15.9))')
127         &       'grdchk output:  ', i, itilepos, jtilepos, layer,
128         &       xxmemo_ref, xxmemo_pert
129            call print_message( msgbuf, standardmessageunit,
130         &                      SQUEEZE_RIGHT , mythid)
131            if ( ierr .eq. 0 ) then
132              write(msgbuf,'(A,5(1x,D15.9))')
133         &      'grdchk output:  ',
134         &          fcref, fcpert, gfd, adxxmemo, ratio
135              call print_message( msgbuf, standardmessageunit,
136         &                        SQUEEZE_RIGHT , mythid)
137            else
138              if ( ierr .eq. -1 ) then
139                write(msgbuf,'(a)')
140         &      ' Component does not exist (zero)'
141              else if ( ierr .eq. -2 ) then
142                write(msgbuf,'(a)')
143         &      ' Component does not exist (negative)'
144              else if ( ierr .eq. -3 ) then
145                write(msgbuf,'(a)')
146         &      ' Component does not exist (too large)'
147              else if ( ierr .eq. -4 ) then
148                write(msgbuf,'(a)')
149         &      ' Component does not exist (land point)'
150              endif
151              call print_message( msgbuf, standardmessageunit,
152         &                        SQUEEZE_RIGHT , mythid)
153            endif
154            write(msgbuf,'(a)')
155         &  ' '
156            call print_message( msgbuf, standardmessageunit,
157         &                      SQUEEZE_RIGHT , mythid)
158          enddo
159    
160    c--   Print final lines.
161          write(msgbuf,'(a)')
162         &' '
163          call print_message( msgbuf, standardmessageunit,
164         &                    SQUEEZE_RIGHT , mythid)
165          write(msgbuf,'(a)')
166         &'// ======================================================='
167          call print_message( msgbuf, standardmessageunit,
168         &                  SQUEEZE_RIGHT , mythid)
169          write(msgbuf,'(a)')
170         &'// Gradient check results  >>> END <<<'
171          call print_message( msgbuf, standardmessageunit,
172         &                    SQUEEZE_RIGHT , mythid)
173          write(msgbuf,'(a)')
174         &'// ======================================================='
175          call print_message( msgbuf, standardmessageunit,
176         &                    SQUEEZE_RIGHT , mythid)
177          write(msgbuf,'(a)')
178         &' '
179          call print_message( msgbuf, standardmessageunit,
180         &                    SQUEEZE_RIGHT , mythid)
181    
182    #endif /* ALLOW_GRADIENT_CHECK */
183    
184          return
185          end
186    

Legend:
Removed from v.1.1  
changed lines
  Added in v.1.2

  ViewVC Help
Powered by ViewVC 1.1.22