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

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

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


Revision 1.3 - (show annotations) (download)
Thu May 30 22:47:26 2002 UTC (22 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint45d_post, checkpoint45c_post
Changes since 1.2: +16 -10 lines
o modifications to gradient check package (Martin Losch)
  - enable centered differences
  - modified format of standard output

1 C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_print.F,v 1.2.4.1 2002/05/30 22:12:32 heimbach Exp $
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 fcpertplus, fcpertminus
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,6(1x,a15))')
98 & 'grdchk output: ', 'I', 'ITILEPOS', 'JTILEPOS',
99 & 'LAYER', 'X(I)', 'X(I)+/-EPS'
100 call print_message( msgbuf, standardmessageunit,
101 & SQUEEZE_RIGHT , mythid)
102 write(msgbuf,'(a,6(1x,a15))')
103 & 'grdchk output: ', 'FC', 'FC1', 'FC2', 'FC1-FC2/(2*EPS)',
104 & 'GRAD(FC)', '1-FDGRD/ADGRD'
105 call print_message( msgbuf, standardmessageunit,
106 & SQUEEZE_RIGHT , mythid)
107
108 c-- Individual checks.
109 if ( ierr_grdchk .eq. 0 ) then
110 numchecks = ichknum - 1
111 else
112 numchecks = maxgrdchecks
113 endif
114
115 do i = 1, numchecks
116 xxmemo_ref = xxmemref (i)
117 xxmemo_pert = xxmempert (i)
118 adxxmemo = adxxmem (i)
119 fcref = fcrmem (i)
120 fcpertplus = fcppmem (i)
121 fcpertminus = fcppmem (i)
122 gfd = gfdmem (i)
123 ratio = ratiomem(i)
124 itile = bimem (i)
125 jtile = bjmem (i)
126 itilepos = ilocmem (i)
127 jtilepos = jlocmem (i)
128 layer = klocmem (i)
129 icomp = icompmem(i)
130 ierr = ierrmem (i)
131
132 write(msgbuf,'(A,4(I16),2(1x,D15.9))')
133 & 'grdchk output: ', i, itilepos, jtilepos, layer,
134 & xxmemo_ref, xxmemo_pert
135 call print_message( msgbuf, standardmessageunit,
136 & SQUEEZE_RIGHT , mythid)
137 if ( ierr .eq. 0 ) then
138 write(msgbuf,'(A,6(1x,D15.9))')
139 & 'grdchk output: ',
140 & fcref, fcpertplus, fcpertminus, gfd, adxxmemo, ratio
141 call print_message( msgbuf, standardmessageunit,
142 & SQUEEZE_RIGHT , mythid)
143 else
144 if ( ierr .eq. -1 ) then
145 write(msgbuf,'(a)')
146 & ' Component does not exist (zero)'
147 else if ( ierr .eq. -2 ) then
148 write(msgbuf,'(a)')
149 & ' Component does not exist (negative)'
150 else if ( ierr .eq. -3 ) then
151 write(msgbuf,'(a)')
152 & ' Component does not exist (too large)'
153 else if ( ierr .eq. -4 ) then
154 write(msgbuf,'(a)')
155 & ' Component does not exist (land point)'
156 endif
157 call print_message( msgbuf, standardmessageunit,
158 & SQUEEZE_RIGHT , mythid)
159 endif
160 write(msgbuf,'(a)')
161 & ' '
162 call print_message( msgbuf, standardmessageunit,
163 & SQUEEZE_RIGHT , mythid)
164 enddo
165
166 c-- Print final lines.
167 write(msgbuf,'(a)')
168 &' '
169 call print_message( msgbuf, standardmessageunit,
170 & SQUEEZE_RIGHT , mythid)
171 write(msgbuf,'(a)')
172 &'// ======================================================='
173 call print_message( msgbuf, standardmessageunit,
174 & SQUEEZE_RIGHT , mythid)
175 write(msgbuf,'(a)')
176 &'// Gradient check results >>> END <<<'
177 call print_message( msgbuf, standardmessageunit,
178 & SQUEEZE_RIGHT , mythid)
179 write(msgbuf,'(a)')
180 &'// ======================================================='
181 call print_message( msgbuf, standardmessageunit,
182 & SQUEEZE_RIGHT , mythid)
183 write(msgbuf,'(a)')
184 &' '
185 call print_message( msgbuf, standardmessageunit,
186 & SQUEEZE_RIGHT , mythid)
187
188 #endif /* ALLOW_GRADIENT_CHECK */
189
190 return
191 end
192

  ViewVC Help
Powered by ViewVC 1.1.22