1 |
heimbach |
1.6 |
C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_print.F,v 1.2.6.1 2002/05/30 19:55:17 heimbach Exp $ |
2 |
heimbach |
1.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 |
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 |
heimbach |
1.2 |
|
53 |
|
|
integer i |
54 |
|
|
integer itile |
55 |
|
|
integer jtile |
56 |
|
|
integer itilepos |
57 |
|
|
integer jtilepos |
58 |
|
|
integer layer |
59 |
|
|
integer icomp |
60 |
|
|
integer ierr |
61 |
|
|
|
62 |
|
|
integer numchecks |
63 |
|
|
|
64 |
|
|
character*(max_len_mbuf) msgbuf |
65 |
|
|
|
66 |
|
|
c == end of interface == |
67 |
|
|
|
68 |
|
|
c-- Print header. |
69 |
|
|
write(msgbuf,'(a)') |
70 |
|
|
&' ' |
71 |
|
|
call print_message( msgbuf, standardmessageunit, |
72 |
|
|
& SQUEEZE_RIGHT , mythid) |
73 |
|
|
write(msgbuf,'(a)') |
74 |
|
|
&'// =======================================================' |
75 |
|
|
call print_message( msgbuf, standardmessageunit, |
76 |
|
|
& SQUEEZE_RIGHT , mythid) |
77 |
|
|
write(msgbuf,'(a)') |
78 |
|
|
&'// Gradient check results >>> START <<<' |
79 |
|
|
call print_message( msgbuf, standardmessageunit, |
80 |
|
|
& SQUEEZE_RIGHT , mythid) |
81 |
|
|
write(msgbuf,'(a)') |
82 |
|
|
&'// =======================================================' |
83 |
|
|
call print_message( msgbuf, standardmessageunit, |
84 |
|
|
& SQUEEZE_RIGHT , mythid) |
85 |
|
|
write(msgbuf,'(a)') |
86 |
|
|
&' ' |
87 |
|
|
call print_message( msgbuf, standardmessageunit, |
88 |
|
|
& SQUEEZE_RIGHT , mythid) |
89 |
|
|
write(msgbuf,'(a)') |
90 |
|
|
&' ' |
91 |
|
|
call print_message( msgbuf, standardmessageunit, |
92 |
|
|
& SQUEEZE_RIGHT , mythid) |
93 |
|
|
|
94 |
|
|
write(msgbuf,'(a,e10.3)') |
95 |
|
|
&' EPS = ',grdchk_eps |
96 |
|
|
call print_message( msgbuf, standardmessageunit, |
97 |
|
|
& SQUEEZE_RIGHT , mythid) |
98 |
|
|
|
99 |
heimbach |
1.6 |
write(msgbuf,'(a,7(1x,a15))') |
100 |
|
|
& 'grdchk output: ', 'procId', 'I', 'ITILEPOS', 'JTILEPOS', |
101 |
heimbach |
1.3 |
& 'LAYER', 'X(I)', 'X(I)+/-EPS' |
102 |
|
|
call print_message( msgbuf, standardmessageunit, |
103 |
|
|
& SQUEEZE_RIGHT , mythid) |
104 |
heimbach |
1.6 |
#ifdef ALLOW_TANGENTLINEAR_RUN |
105 |
|
|
write(msgbuf,'(a,7(1x,a15))') |
106 |
|
|
& 'grdchk output: ', ' ', 'FC', 'FC1', 'FC2', |
107 |
|
|
& 'FC1-FC2/(2*EPS)', 'TLM GRAD(FC)', '1-FDGRD/TLMGRD' |
108 |
|
|
#else |
109 |
|
|
write(msgbuf,'(a,7(1x,a15))') |
110 |
|
|
& 'grdchk output: ', ' ', 'FC', 'FC1', 'FC2', |
111 |
|
|
& 'FC1-FC2/(2*EPS)', 'ADJ GRAD(FC)', '1-FDGRD/ADGRD' |
112 |
|
|
#endif |
113 |
heimbach |
1.2 |
call print_message( msgbuf, standardmessageunit, |
114 |
|
|
& SQUEEZE_RIGHT , mythid) |
115 |
|
|
|
116 |
|
|
c-- Individual checks. |
117 |
|
|
if ( ierr_grdchk .eq. 0 ) then |
118 |
heimbach |
1.6 |
numchecks = ichknum |
119 |
heimbach |
1.2 |
else |
120 |
|
|
numchecks = maxgrdchecks |
121 |
|
|
endif |
122 |
|
|
|
123 |
|
|
do i = 1, numchecks |
124 |
|
|
xxmemo_ref = xxmemref (i) |
125 |
|
|
xxmemo_pert = xxmempert (i) |
126 |
|
|
adxxmemo = adxxmem (i) |
127 |
heimbach |
1.5 |
ftlxxmemo = ftlxxmem (i) |
128 |
heimbach |
1.3 |
fcref = fcrmem (i) |
129 |
|
|
fcpertplus = fcppmem (i) |
130 |
heimbach |
1.6 |
fcpertminus = fcpmmem (i) |
131 |
|
|
gfd = gfdmem (i) |
132 |
heimbach |
1.5 |
ratio_ad = ratioadmem (i) |
133 |
|
|
ratio_ftl = ratioftlmem (i) |
134 |
|
|
itile = bimem (i) |
135 |
|
|
jtile = bjmem (i) |
136 |
|
|
itilepos = ilocmem (i) |
137 |
|
|
jtilepos = jlocmem (i) |
138 |
|
|
layer = klocmem (i) |
139 |
|
|
icomp = icompmem(i) |
140 |
|
|
ierr = ierrmem (i) |
141 |
heimbach |
1.2 |
|
142 |
heimbach |
1.6 |
write(msgbuf,'(A,5(I16),2(1x,D15.9))') |
143 |
|
|
& 'grdchk output: ', |
144 |
|
|
& myprocid, i, itilepos, jtilepos, layer, |
145 |
heimbach |
1.2 |
& xxmemo_ref, xxmemo_pert |
146 |
|
|
call print_message( msgbuf, standardmessageunit, |
147 |
|
|
& SQUEEZE_RIGHT , mythid) |
148 |
|
|
if ( ierr .eq. 0 ) then |
149 |
heimbach |
1.6 |
#ifdef ALLOW_TANGENTLINEAR_RUN |
150 |
|
|
write(msgbuf,'(A,1(1x,a15),6(1x,D15.9))') |
151 |
|
|
& 'grdchk output: ', ' ', |
152 |
|
|
& fcref, fcpertplus, fcpertminus, |
153 |
|
|
& gfd, ftlxxmemo, ratio_ftl |
154 |
|
|
#else |
155 |
|
|
write(msgbuf,'(A,1(1x,a15),6(1x,D15.9))') |
156 |
|
|
& 'grdchk output: ', ' ', |
157 |
|
|
& fcref, fcpertplus, fcpertminus, |
158 |
|
|
& gfd, adxxmemo, ratio_ad |
159 |
|
|
#endif |
160 |
heimbach |
1.2 |
call print_message( msgbuf, standardmessageunit, |
161 |
|
|
& SQUEEZE_RIGHT , mythid) |
162 |
|
|
else |
163 |
|
|
if ( ierr .eq. -1 ) then |
164 |
|
|
write(msgbuf,'(a)') |
165 |
|
|
& ' Component does not exist (zero)' |
166 |
|
|
else if ( ierr .eq. -2 ) then |
167 |
|
|
write(msgbuf,'(a)') |
168 |
|
|
& ' Component does not exist (negative)' |
169 |
|
|
else if ( ierr .eq. -3 ) then |
170 |
|
|
write(msgbuf,'(a)') |
171 |
|
|
& ' Component does not exist (too large)' |
172 |
|
|
else if ( ierr .eq. -4 ) then |
173 |
|
|
write(msgbuf,'(a)') |
174 |
|
|
& ' Component does not exist (land point)' |
175 |
|
|
endif |
176 |
|
|
call print_message( msgbuf, standardmessageunit, |
177 |
|
|
& SQUEEZE_RIGHT , mythid) |
178 |
|
|
endif |
179 |
|
|
write(msgbuf,'(a)') |
180 |
|
|
& ' ' |
181 |
|
|
call print_message( msgbuf, standardmessageunit, |
182 |
|
|
& SQUEEZE_RIGHT , mythid) |
183 |
|
|
enddo |
184 |
|
|
|
185 |
|
|
c-- Print final lines. |
186 |
|
|
write(msgbuf,'(a)') |
187 |
|
|
&' ' |
188 |
|
|
call print_message( msgbuf, standardmessageunit, |
189 |
|
|
& SQUEEZE_RIGHT , mythid) |
190 |
|
|
write(msgbuf,'(a)') |
191 |
|
|
&'// =======================================================' |
192 |
|
|
call print_message( msgbuf, standardmessageunit, |
193 |
|
|
& SQUEEZE_RIGHT , mythid) |
194 |
|
|
write(msgbuf,'(a)') |
195 |
|
|
&'// Gradient check results >>> END <<<' |
196 |
|
|
call print_message( msgbuf, standardmessageunit, |
197 |
|
|
& SQUEEZE_RIGHT , mythid) |
198 |
|
|
write(msgbuf,'(a)') |
199 |
|
|
&'// =======================================================' |
200 |
|
|
call print_message( msgbuf, standardmessageunit, |
201 |
|
|
& SQUEEZE_RIGHT , mythid) |
202 |
|
|
write(msgbuf,'(a)') |
203 |
|
|
&' ' |
204 |
|
|
call print_message( msgbuf, standardmessageunit, |
205 |
|
|
& SQUEEZE_RIGHT , mythid) |
206 |
|
|
|
207 |
|
|
#endif /* ALLOW_GRADIENT_CHECK */ |
208 |
|
|
|
209 |
|
|
return |
210 |
|
|
end |
211 |
|
|
|