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

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

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


Revision 1.13 - (hide annotations) (download)
Fri Aug 5 17:56:55 2011 UTC (12 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63g, checkpoint63h, checkpoint63i, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63b, checkpoint63c
Changes since 1.12: +6 -6 lines
also print (local) tile indices

1 jmc 1.13 C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_print.F,v 1.12 2011/05/24 22:40:30 jmc Exp $
2 heimbach 1.9 C $Name: $
3 heimbach 1.2
4 jmc 1.12 #include "GRDCHK_OPTIONS.h"
5 heimbach 1.2
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 jmc 1.10 #include "SIZE.h"
30 heimbach 1.2 #include "EEPARAMS.h"
31     #include "grdchk.h"
32    
33     c == routine arguments ==
34    
35     integer ichknum
36     integer ierr_grdchk
37 jmc 1.10 integer mythid
38 heimbach 1.2
39 heimbach 1.9 #ifdef ALLOW_GRDCHK
40 heimbach 1.2 c == local variables ==
41    
42     _RL fcref
43 heimbach 1.3 _RL fcpertplus, fcpertminus
44 heimbach 1.2 _RL xxmemo_ref
45     _RL xxmemo_pert
46     _RL gfd
47     _RL adxxmemo
48 heimbach 1.5 _RL ftlxxmemo
49     _RL ratio_ad
50     _RL ratio_ftl
51 heimbach 1.2
52     integer i
53     integer itile
54     integer jtile
55     integer itilepos
56     integer jtilepos
57     integer layer
58     integer icomp
59     integer ierr
60    
61     integer numchecks
62    
63     character*(max_len_mbuf) msgbuf
64    
65     c == end of interface ==
66    
67     c-- Print header.
68     write(msgbuf,'(a)')
69     &' '
70     call print_message( msgbuf, standardmessageunit,
71 jmc 1.10 & SQUEEZE_RIGHT, mythid )
72 heimbach 1.2 write(msgbuf,'(a)')
73     &'// ======================================================='
74     call print_message( msgbuf, standardmessageunit,
75 jmc 1.10 & SQUEEZE_RIGHT, mythid )
76 heimbach 1.2 write(msgbuf,'(a)')
77     &'// Gradient check results >>> START <<<'
78     call print_message( msgbuf, standardmessageunit,
79 jmc 1.10 & SQUEEZE_RIGHT, mythid)
80 heimbach 1.2 write(msgbuf,'(a)')
81     &'// ======================================================='
82     call print_message( msgbuf, standardmessageunit,
83 jmc 1.10 & SQUEEZE_RIGHT , mythid )
84 heimbach 1.2 write(msgbuf,'(a)')
85     &' '
86     call print_message( msgbuf, standardmessageunit,
87 jmc 1.10 & SQUEEZE_RIGHT, mythid )
88    
89     write(msgbuf,'(A,1PE14.6)')
90     &' EPS = ',grdchk_eps
91     call print_message( msgbuf, standardmessageunit,
92     & SQUEEZE_RIGHT, mythid )
93 heimbach 1.2 write(msgbuf,'(a)')
94     &' '
95     call print_message( msgbuf, standardmessageunit,
96 jmc 1.10 & SQUEEZE_RIGHT, mythid )
97 heimbach 1.2
98 jmc 1.13 write(msgbuf,'(A,2X,4A,3(3X,A),11X,A)')
99     & 'grdchk output h.p:', 'Id', ' Itile', ' Jtile',
100     & ' LAYER', 'bi', 'bj', 'X(Id)', 'X(Id)+/-EPS'
101 heimbach 1.2 call print_message( msgbuf, standardmessageunit,
102 jmc 1.10 & SQUEEZE_RIGHT , mythid )
103 jmc 1.11 write(msgbuf,'(A,2X,A,A4,1X,2A21)')
104 jmc 1.10 & 'grdchk output h.c:', 'Id', 'FC', 'FC1', 'FC2'
105 heimbach 1.3 call print_message( msgbuf, standardmessageunit,
106 jmc 1.10 & SQUEEZE_RIGHT, mythid )
107 heimbach 1.6 #ifdef ALLOW_TANGENTLINEAR_RUN
108 jmc 1.11 write(msgbuf,'(A,2X,A,2X,2A18,4X,A18)')
109 jmc 1.10 & 'grdchk output h.g:', 'Id',
110 heimbach 1.6 & 'FC1-FC2/(2*EPS)', 'TLM GRAD(FC)', '1-FDGRD/TLMGRD'
111     #else
112 jmc 1.11 write(msgbuf,'(A,2X,A,2X,2A18,4X,A18)')
113 jmc 1.10 & 'grdchk output h.g:', 'Id',
114 heimbach 1.6 & 'FC1-FC2/(2*EPS)', 'ADJ GRAD(FC)', '1-FDGRD/ADGRD'
115     #endif
116 heimbach 1.2 call print_message( msgbuf, standardmessageunit,
117 jmc 1.10 & SQUEEZE_RIGHT, mythid )
118 heimbach 1.2
119     c-- Individual checks.
120     if ( ierr_grdchk .eq. 0 ) then
121 heimbach 1.6 numchecks = ichknum
122 heimbach 1.2 else
123     numchecks = maxgrdchecks
124     endif
125    
126     do i = 1, numchecks
127     xxmemo_ref = xxmemref (i)
128     xxmemo_pert = xxmempert (i)
129     adxxmemo = adxxmem (i)
130 heimbach 1.5 ftlxxmemo = ftlxxmem (i)
131 heimbach 1.3 fcref = fcrmem (i)
132     fcpertplus = fcppmem (i)
133 heimbach 1.6 fcpertminus = fcpmmem (i)
134     gfd = gfdmem (i)
135 heimbach 1.5 ratio_ad = ratioadmem (i)
136     ratio_ftl = ratioftlmem (i)
137     itile = bimem (i)
138     jtile = bjmem (i)
139     itilepos = ilocmem (i)
140     jtilepos = jlocmem (i)
141     layer = klocmem (i)
142     icomp = icompmem(i)
143     ierr = ierrmem (i)
144 heimbach 1.2
145 jmc 1.10 write(msgbuf,'(a)')
146     & ' '
147     call print_message( msgbuf, standardmessageunit,
148     & SQUEEZE_RIGHT, mythid )
149 jmc 1.13 write(msgbuf,'(A,I4,3I6,2I5,1x,1P2E17.9)')
150 jmc 1.10 & 'grdchk output (p):',
151 jmc 1.13 & i, itilepos, jtilepos, layer, itile, jtile,
152 heimbach 1.2 & xxmemo_ref, xxmemo_pert
153     call print_message( msgbuf, standardmessageunit,
154 jmc 1.10 & SQUEEZE_RIGHT, mythid )
155 heimbach 1.2 if ( ierr .eq. 0 ) then
156 jmc 1.11 write(msgbuf,'(A,I4,1P3E21.13)')
157 jmc 1.10 & 'grdchk output (c):',
158     & i, fcref, fcpertplus, fcpertminus
159     call print_message( msgbuf, standardmessageunit,
160     & SQUEEZE_RIGHT, mythid )
161 heimbach 1.6 #ifdef ALLOW_TANGENTLINEAR_RUN
162 jmc 1.11 write(msgbuf,'(A,I4,3x,1P3E21.13)')
163 jmc 1.10 & 'grdchk output (g):',
164     & i, gfd, ftlxxmemo, ratio_ftl
165 heimbach 1.6 #else
166 jmc 1.11 write(msgbuf,'(A,I4,3x,1P3E21.13)')
167 jmc 1.10 & 'grdchk output (g):',
168     & i, gfd, adxxmemo, ratio_ad
169 heimbach 1.6 #endif
170 heimbach 1.2 call print_message( msgbuf, standardmessageunit,
171 jmc 1.10 & SQUEEZE_RIGHT, mythid )
172 heimbach 1.2 else
173     if ( ierr .eq. -1 ) then
174     write(msgbuf,'(a)')
175     & ' Component does not exist (zero)'
176     else if ( ierr .eq. -2 ) then
177     write(msgbuf,'(a)')
178     & ' Component does not exist (negative)'
179     else if ( ierr .eq. -3 ) then
180     write(msgbuf,'(a)')
181     & ' Component does not exist (too large)'
182     else if ( ierr .eq. -4 ) then
183     write(msgbuf,'(a)')
184     & ' Component does not exist (land point)'
185     endif
186     call print_message( msgbuf, standardmessageunit,
187 jmc 1.10 & SQUEEZE_RIGHT, mythid )
188 heimbach 1.2 endif
189     enddo
190    
191     c-- Print final lines.
192     write(msgbuf,'(a)')
193     &' '
194     call print_message( msgbuf, standardmessageunit,
195 jmc 1.10 & SQUEEZE_RIGHT, mythid )
196 heimbach 1.2 write(msgbuf,'(a)')
197     &'// ======================================================='
198     call print_message( msgbuf, standardmessageunit,
199 jmc 1.10 & SQUEEZE_RIGHT, mythid )
200 heimbach 1.2 write(msgbuf,'(a)')
201     &'// Gradient check results >>> END <<<'
202     call print_message( msgbuf, standardmessageunit,
203 jmc 1.10 & SQUEEZE_RIGHT, mythid )
204 heimbach 1.2 write(msgbuf,'(a)')
205     &'// ======================================================='
206     call print_message( msgbuf, standardmessageunit,
207 jmc 1.10 & SQUEEZE_RIGHT, mythid )
208 heimbach 1.2 write(msgbuf,'(a)')
209     &' '
210     call print_message( msgbuf, standardmessageunit,
211 jmc 1.10 & SQUEEZE_RIGHT, mythid )
212 heimbach 1.2
213 heimbach 1.9 #endif /* ALLOW_GRDCHK */
214 heimbach 1.2
215     return
216     end

  ViewVC Help
Powered by ViewVC 1.1.22