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