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

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

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


Revision 1.40 - (hide annotations) (download)
Wed Aug 22 19:34:51 2012 UTC (11 years, 9 months ago) by heimbach
Branch: MAIN
Changes since 1.39: +20 -4 lines
Merge OpenAD specific stuff to pkg/grdchk

1 heimbach 1.40 C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_main.F,v 1.39 2012/08/21 19:49:41 gforget Exp $
2 heimbach 1.11 C $Name: $
3 heimbach 1.2
4 jmc 1.31 #include "GRDCHK_OPTIONS.h"
5 jmc 1.36 #include "AD_CONFIG.h"
6 heimbach 1.2
7 heimbach 1.3 CBOI
8     C
9     C !TITLE: GRADIENT CHECK
10     C !AUTHORS: mitgcm developers ( support@mitgcm.org )
11     C !AFFILIATION: Massachussetts Institute of Technology
12     C !DATE:
13     C !INTRODUCTION: gradient check package
14     c \bv
15     c Compare the gradients calculated by the adjoint model with
16     c finite difference approximations.
17     c
18     C !CALLING SEQUENCE:
19     c
20     c the_model_main
21     c |
22     c |-- ctrl_unpack
23     c |-- adthe_main_loop - unperturbed cost function and
24     c |-- ctrl_pack adjoint gradient are computed here
25     c |
26     c |-- grdchk_main
27     c |
28     c |-- grdchk_init
29     c |-- do icomp=... - loop over control vector elements
30     c |
31     c |-- grdchk_loc - determine location of icomp on grid
32     c |
33     c |-- grdchk_getxx - get control vector component from file
34     c | perturb it and write back to file
35 jmc 1.31 c |-- grdchk_getadxx - get gradient component calculated
36 heimbach 1.3 c | via adjoint
37     c |-- the_main_loop - forward run and cost evaluation
38     c | with perturbed control vector element
39     c |-- calculate ratio of adj. vs. finite difference gradient
40     c |
41     c |-- grdchk_setxx - Reset control vector element
42     c |
43     c |-- grdchk_print - print results
44     c \ev
45     CEOI
46    
47     CBOP
48     C !ROUTINE: grdchk_main
49     C !INTERFACE:
50     subroutine grdchk_main( mythid )
51 heimbach 1.2
52 heimbach 1.3 C !DESCRIPTION: \bv
53 heimbach 1.2 c ==================================================================
54     c SUBROUTINE grdchk_main
55     c ==================================================================
56     c o Compare the gradients calculated by the adjoint model with
57     c finite difference approximations.
58     c started: Christian Eckert eckert@mit.edu 24-Feb-2000
59 heimbach 1.4 c continued&finished: heimbach@mit.edu: 13-Jun-2001
60     c changed: mlosch@ocean.mit.edu: 09-May-2002
61     c - added centered difference vs. 1-sided difference option
62     c - improved output format for readability
63     c - added control variable hFacC
64 heimbach 1.7 c heimbach@mit.edu 24-Feb-2003
65     c - added tangent linear gradient checks
66     c - fixes for multiproc. gradient checks
67     c - added more control variables
68 jmc 1.31 c
69 heimbach 1.2 c ==================================================================
70     c SUBROUTINE grdchk_main
71     c ==================================================================
72 heimbach 1.3 C \ev
73 heimbach 1.2
74 heimbach 1.3 C !USES:
75 heimbach 1.2 implicit none
76    
77     c == global variables ==
78     #include "SIZE.h"
79     #include "EEPARAMS.h"
80     #include "PARAMS.h"
81     #include "grdchk.h"
82     #include "cost.h"
83 heimbach 1.21 #include "ctrl.h"
84 heimbach 1.9 #ifdef ALLOW_TANGENTLINEAR_RUN
85     #include "g_cost.h"
86     #endif
87 heimbach 1.2
88 heimbach 1.3 C !INPUT/OUTPUT PARAMETERS:
89 heimbach 1.2 c == routine arguments ==
90 jmc 1.31 integer mythid
91 heimbach 1.2
92 heimbach 1.11 #ifdef ALLOW_GRDCHK
93 heimbach 1.3 C !LOCAL VARIABLES:
94 heimbach 1.2 c == local variables ==
95     integer myiter
96 jmc 1.31 _RL mytime
97 heimbach 1.2 integer bi, itlo, ithi
98     integer bj, jtlo, jthi
99     integer i, imin, imax
100     integer j, jmin, jmax
101     integer k
102    
103     integer icomp
104     integer ichknum
105     integer icvrec
106     integer jtile
107     integer itile
108     integer layer
109 heimbach 1.8 integer obcspos
110 heimbach 1.2 integer itilepos
111     integer jtilepos
112 heimbach 1.19 integer icglo
113 heimbach 1.2 integer itest
114     integer ierr
115     integer ierr_grdchk
116     _RL gfd
117     _RL fcref
118 heimbach 1.4 _RL fcpertplus, fcpertminus
119 heimbach 1.6 _RL ratio_ad
120     _RL ratio_ftl
121 heimbach 1.2 _RL xxmemo_ref
122     _RL xxmemo_pert
123     _RL adxxmemo
124 heimbach 1.6 _RL ftlxxmemo
125     _RL localEps
126     _RL grdchk_epsfac
127    
128 heimbach 1.7 _RL tmpplot1(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy)
129     _RL tmpplot2(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy)
130     _RL tmpplot3(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy)
131    
132 heimbach 1.12 CHARACTER*(MAX_LEN_MBUF) msgBuf
133    
134 heimbach 1.2 c == end of interface ==
135 heimbach 1.3 CEOP
136 heimbach 1.2
137     c-- Set the loop ranges.
138     jtlo = mybylo(mythid)
139     jthi = mybyhi(mythid)
140     itlo = mybxlo(mythid)
141     ithi = mybxhi(mythid)
142     jmin = 1
143     jmax = sny
144     imin = 1
145     imax = snx
146    
147 heimbach 1.16 print *, 'ph-check entering grdchk_main '
148    
149 heimbach 1.2 c-- initialise variables
150     call grdchk_init( mythid )
151    
152 jmc 1.29 c-- Compute the adjoint model gradients.
153 heimbach 1.2 c-- Compute the unperturbed cost function.
154 heimbach 1.7 cph Gradient via adjoint has already been computed,
155     cph and so has unperturbed cost function,
156     cph assuming all xx_ fields are initialised to zero.
157 heimbach 1.2
158 heimbach 1.28 ierr = 0
159 heimbach 1.9 ierr_grdchk = 0
160 heimbach 1.27 adxxmemo = 0.
161     ftlxxmemo = 0.
162 heimbach 1.40 #if (defined (ALLOW_ADMTLM))
163 heimbach 1.21 fcref = objf_state_final(idep,jdep,1,1,1)
164 heimbach 1.40 #elif (defined (ALLOW_AUTODIFF_OPENAD))
165     fcref = fc%v
166 heimbach 1.21 #else
167 heimbach 1.2 fcref = fc
168 heimbach 1.21 #endif
169 heimbach 1.9
170     print *, 'ph-check fcref = ', fcref
171 heimbach 1.2
172     do bj = jtlo, jthi
173     do bi = itlo, ithi
174     do k = 1, nr
175     do j = jmin, jmax
176     do i = imin, imax
177     tmpplot1(i,j,k,bi,bj) = 0. _d 0
178     tmpplot2(i,j,k,bi,bj) = 0. _d 0
179 heimbach 1.6 tmpplot3(i,j,k,bi,bj) = 0. _d 0
180 heimbach 1.2 end do
181     end do
182     end do
183     end do
184     end do
185    
186 heimbach 1.4 if ( useCentralDiff ) then
187     grdchk_epsfac = 2. _d 0
188     else
189     grdchk_epsfac = 1. _d 0
190     end if
191    
192 heimbach 1.23 WRITE(standardmessageunit,'(A)')
193     & 'grad-res -------------------------------'
194     WRITE(standardmessageunit,'(2a)')
195 heimbach 1.19 & ' grad-res proc # i j k bi bj iobc',
196 heimbach 1.18 & ' fc ref fc + eps fc - eps'
197 heimbach 1.7 #ifdef ALLOW_TANGENTLINEAR_RUN
198 heimbach 1.23 WRITE(standardmessageunit,'(2a)')
199 heimbach 1.19 & ' grad-res proc # i j k bi bj iobc',
200 heimbach 1.18 & ' tlm grad fd grad 1 - fd/tlm'
201 heimbach 1.7 #else
202 heimbach 1.23 WRITE(standardmessageunit,'(2a)')
203 heimbach 1.19 & ' grad-res proc # i j k bi bj iobc',
204 heimbach 1.18 & ' adj grad fd grad 1 - fd/adj'
205 heimbach 1.7 #endif
206    
207 heimbach 1.2 c-- Compute the finite difference approximations.
208     c-- Cycle through all processes doing NINT(nend-nbeg+1)/nstep
209     c-- gradient checks.
210 heimbach 1.14
211 jmc 1.31 if ( nbeg .EQ. 0 )
212 heimbach 1.17 & call grdchk_get_position( mythid )
213 heimbach 1.2
214 heimbach 1.7 do icomp = nbeg, nend, nstep
215 heimbach 1.2
216 heimbach 1.7 ichknum = (icomp - nbeg)/nstep + 1
217 jmc 1.32 adxxmemo = 0.
218 heimbach 1.2
219 heimbach 1.19 cph(
220 jmc 1.31 cph-print print *, 'ph-grd _main: nbeg, icomp, ichknum ',
221 heimbach 1.19 cph-print & nbeg, icomp, ichknum
222     cph)
223 heimbach 1.7 if (ichknum .le. maxgrdchecks ) then
224 heimbach 1.2
225 heimbach 1.7 c-- Determine the location of icomp on the grid.
226     if ( myProcId .EQ. grdchkwhichproc ) then
227 jmc 1.31 call grdchk_loc( icomp, ichknum,
228 heimbach 1.8 & icvrec, itile, jtile, layer, obcspos,
229 heimbach 1.19 & itilepos, jtilepos, icglo, itest, ierr,
230 heimbach 1.7 & mythid )
231 heimbach 1.18 cph(
232 heimbach 1.19 cph-print print *, 'ph-grd ----- back from loc -----',
233     cph-print & icvrec, itilepos, jtilepos, layer, obcspos
234 heimbach 1.18 cph)
235 jmc 1.34 else
236 jmc 1.35 icvrec = 0
237 jmc 1.34 itile = 0
238     jtile = 0
239     layer = 0
240 jmc 1.35 obcspos = 0
241 jmc 1.34 itilepos = 0
242     jtilepos = 0
243 jmc 1.35 icglo = 0
244     itest = 0
245 heimbach 1.7 endif
246 jmc 1.31
247 gforget 1.38 c make sure that all procs have correct file records, so that useSingleCpuIO works
248     CALL GLOBAL_SUM_INT( icvrec , myThid )
249     CALL GLOBAL_SUM_INT( layer , myThid )
250 gforget 1.39 CALL GLOBAL_SUM_INT( ierr , myThid )
251 gforget 1.38
252 heimbach 1.6 c******************************************************
253     c-- (A): get gradient component calculated via adjoint
254     c******************************************************
255 heimbach 1.7
256     c-- get gradient component calculated via adjoint
257     call grdchk_getadxx( icvrec,
258     & itile, jtile, layer,
259     & itilepos, jtilepos,
260 gforget 1.38 & adxxmemo, ierr, mythid )
261 jmc 1.32 C-- Add a global-sum call so that all proc will get the adjoint gradient
262     _GLOBAL_SUM_RL( adxxmemo, myThid )
263 heimbach 1.4
264 heimbach 1.6 #ifdef ALLOW_TANGENTLINEAR_RUN
265     c******************************************************
266     c-- (B): Get gradient component g_fc from tangent linear run:
267     c******************************************************
268     c--
269     c-- 1. perturb control vector component: xx(i)=1.
270    
271 heimbach 1.7 localEps = 1. _d 0
272     call grdchk_getxx( icvrec, TANGENT_SIMULATION,
273     & itile, jtile, layer,
274     & itilepos, jtilepos,
275     & xxmemo_ref, xxmemo_pert, localEps,
276 gforget 1.38 & ierr, mythid )
277 heimbach 1.6
278     c--
279     c-- 2. perform tangent linear run
280 heimbach 1.7 mytime = starttime
281     myiter = niter0
282 heimbach 1.21 #ifdef ALLOW_ADMTLM
283     do k=1,4*Nr+1
284     do j=1,sny
285     do i=1,snx
286     g_objf_state_final(i,j,1,1,k) = 0.
287     enddo
288     enddo
289     enddo
290     #else
291 heimbach 1.7 g_fc = 0.
292 heimbach 1.21 #endif
293    
294 heimbach 1.7 call g_the_main_loop( mytime, myiter, mythid )
295 heimbach 1.21 #ifdef ALLOW_ADMTLM
296     ftlxxmemo = g_objf_state_final(idep,jdep,1,1,1)
297     #else
298 heimbach 1.7 ftlxxmemo = g_fc
299 heimbach 1.21 #endif
300    
301 heimbach 1.6 c--
302     c-- 3. reset control vector
303 heimbach 1.7 call grdchk_setxx( icvrec, TANGENT_SIMULATION,
304     & itile, jtile, layer,
305     & itilepos, jtilepos,
306 gforget 1.38 & xxmemo_ref, ierr, mythid )
307 heimbach 1.7
308     #endif /* ALLOW_TANGENTLINEAR_RUN */
309    
310 heimbach 1.6
311     c******************************************************
312     c-- (C): Get gradient via finite difference perturbation
313     c******************************************************
314    
315 heimbach 1.2 c-- get control vector component from file
316 heimbach 1.7 c-- perturb it and write back to file
317 heimbach 1.6 c-- positive perturbation
318 heimbach 1.7 localEps = abs(grdchk_eps)
319 gforget 1.38 call grdchk_getxx( icvrec, FORWARD_SIMULATION,
320 heimbach 1.7 & itile, jtile, layer,
321     & itilepos, jtilepos,
322     & xxmemo_ref, xxmemo_pert, localEps,
323 gforget 1.38 & ierr, mythid )
324 jmc 1.31
325 heimbach 1.4 c-- forward run with perturbed control vector
326 heimbach 1.7 mytime = starttime
327     myiter = niter0
328 heimbach 1.40 #ifdef ALLOW_AUTODIFF_OPENAD
329     call OpenAD_the_main_loop( mytime, myiter, mythid )
330     #else
331 heimbach 1.7 call the_main_loop( mytime, myiter, mythid )
332 heimbach 1.40 #endif
333    
334     #if (defined (ALLOW_ADMTLM))
335 heimbach 1.21 fcpertplus = objf_state_final(idep,jdep,1,1,1)
336 heimbach 1.40 #elif (defined (ALLOW_AUTODIFF_OPENAD))
337     fcpertplus = fc%v
338 heimbach 1.21 #else
339 heimbach 1.7 fcpertplus = fc
340 heimbach 1.21 #endif
341 jmc 1.37 print *, 'ph-check fcpertplus = ', fcpertplus
342 jmc 1.31
343 heimbach 1.4 c-- Reset control vector.
344 gforget 1.38 call grdchk_setxx( icvrec, FORWARD_SIMULATION,
345 heimbach 1.7 & itile, jtile, layer,
346     & itilepos, jtilepos,
347 gforget 1.38 & xxmemo_ref, ierr, mythid )
348 heimbach 1.2
349 heimbach 1.7 fcpertminus = fcref
350 heimbach 1.13 print *, 'ph-check fcpertminus = ', fcpertminus
351 heimbach 1.4
352 heimbach 1.7 if ( useCentralDiff ) then
353 heimbach 1.4
354 heimbach 1.6 c-- get control vector component from file
355 heimbach 1.7 c-- perturb it and write back to file
356 heimbach 1.4 c-- repeat the proceedure for a negative perturbation
357 gforget 1.38 localEps = - abs(grdchk_eps)
358     call grdchk_getxx( icvrec, FORWARD_SIMULATION,
359 heimbach 1.7 & itile, jtile, layer,
360     & itilepos, jtilepos,
361     & xxmemo_ref, xxmemo_pert, localEps,
362 gforget 1.38 & ierr, mythid )
363 jmc 1.31
364 heimbach 1.2 c-- forward run with perturbed control vector
365 heimbach 1.7 mytime = starttime
366     myiter = niter0
367 heimbach 1.40 #ifdef ALLOW_AUTODIFF_OPENAD
368     call OpenAD_the_main_loop( mytime, myiter, mythid )
369     #else
370 heimbach 1.7 call the_main_loop( mytime, myiter, mythid )
371 heimbach 1.40 #endif
372    
373     #if (defined (ALLOW_ADMTLM))
374 heimbach 1.21 fcpertminus = objf_state_final(idep,jdep,1,1,1)
375 heimbach 1.40 #elif (defined (ALLOW_AUTODIFF_OPENAD))
376     fcpertminus = fc%v
377 heimbach 1.21 #else
378 heimbach 1.7 fcpertminus = fc
379 heimbach 1.21 #endif
380 jmc 1.31
381 heimbach 1.4 c-- Reset control vector.
382 gforget 1.38 call grdchk_setxx( icvrec, FORWARD_SIMULATION,
383 heimbach 1.7 & itile, jtile, layer,
384     & itilepos, jtilepos,
385 gforget 1.38 & xxmemo_ref, ierr, mythid )
386 heimbach 1.7
387     c-- end of if useCentralDiff ...
388     end if
389 heimbach 1.4
390 heimbach 1.6 c******************************************************
391     c-- (D): calculate relative differences between gradients
392     c******************************************************
393    
394 jmc 1.33 if ( grdchk_eps .eq. 0. ) then
395     gfd = (fcpertplus-fcpertminus)
396     else
397     gfd = (fcpertplus-fcpertminus)
398     & /(grdchk_epsfac*grdchk_eps)
399     endif
400 heimbach 1.2
401 jmc 1.33 if ( adxxmemo .eq. 0. ) then
402     ratio_ad = abs( adxxmemo - gfd )
403     else
404     ratio_ad = 1. - gfd/adxxmemo
405     endif
406 jmc 1.31
407 jmc 1.33 if ( ftlxxmemo .eq. 0. ) then
408     ratio_ftl = abs( ftlxxmemo - gfd )
409     else
410     ratio_ftl = 1. - gfd/ftlxxmemo
411     endif
412 heimbach 1.7
413 jmc 1.33 if ( myProcId .EQ. grdchkwhichproc .AND.
414     & ierr .EQ. 0 ) then
415 heimbach 1.7 tmpplot1(itilepos,jtilepos,layer,itile,jtile)
416     & = gfd
417     tmpplot2(itilepos,jtilepos,layer,itile,jtile)
418     & = ratio_ad
419     tmpplot3(itilepos,jtilepos,layer,itile,jtile)
420     & = ratio_ftl
421 jmc 1.34 endif
422 heimbach 1.7
423 jmc 1.34 if ( ierr .EQ. 0 ) then
424 heimbach 1.7 fcrmem ( ichknum ) = fcref
425     fcppmem ( ichknum ) = fcpertplus
426     fcpmmem ( ichknum ) = fcpertminus
427     xxmemref ( ichknum ) = xxmemo_ref
428     xxmempert ( ichknum ) = xxmemo_pert
429     gfdmem ( ichknum ) = gfd
430     adxxmem ( ichknum ) = adxxmemo
431     ftlxxmem ( ichknum ) = ftlxxmemo
432     ratioadmem ( ichknum ) = ratio_ad
433     ratioftlmem ( ichknum ) = ratio_ftl
434    
435     irecmem ( ichknum ) = icvrec
436     bimem ( ichknum ) = itile
437     bjmem ( ichknum ) = jtile
438     ilocmem ( ichknum ) = itilepos
439     jlocmem ( ichknum ) = jtilepos
440     klocmem ( ichknum ) = layer
441 heimbach 1.8 iobcsmem ( ichknum ) = obcspos
442 heimbach 1.7 icompmem ( ichknum ) = icomp
443     ichkmem ( ichknum ) = ichknum
444     itestmem ( ichknum ) = itest
445     ierrmem ( ichknum ) = ierr
446 heimbach 1.19 icglomem ( ichknum ) = icglo
447 jmc 1.34 endif
448    
449     if ( myProcId .EQ. grdchkwhichproc .AND.
450     & ierr .EQ. 0 ) then
451 heimbach 1.19
452 heimbach 1.23 WRITE(standardmessageunit,'(A)')
453     & 'grad-res -------------------------------'
454 jmc 1.32 WRITE(standardmessageunit,'(A,8I5,1x,1P3E19.11)')
455 heimbach 1.23 & ' grad-res ',myprocid,ichknum,itilepos,jtilepos,
456 heimbach 1.19 & layer,itile,jtile,obcspos,
457 heimbach 1.7 & fcref, fcpertplus, fcpertminus
458     #ifdef ALLOW_TANGENTLINEAR_RUN
459 jmc 1.32 WRITE(standardmessageunit,'(A,8I5,1x,1P3E19.11)')
460 heimbach 1.23 & ' grad-res ',myprocid,ichknum,ichkmem(ichknum),
461 heimbach 1.19 & icompmem(ichknum),itestmem(ichknum),
462     & bimem(ichknum),bjmem(ichknum),iobcsmem(ichknum),
463 heimbach 1.7 & ftlxxmemo, gfd, ratio_ftl
464     #else
465 jmc 1.32 WRITE(standardmessageunit,'(A,8I5,1x,1P3E19.11)')
466 heimbach 1.23 & ' grad-res ',myprocid,ichknum,ichkmem(ichknum),
467 heimbach 1.19 & icompmem(ichknum),itestmem(ichknum),
468     & bimem(ichknum),bjmem(ichknum),obcspos,
469 heimbach 1.7 & adxxmemo, gfd, ratio_ad
470 jmc 1.32 #endif
471     endif
472     #ifdef ALLOW_TANGENTLINEAR_RUN
473 jmc 1.37 WRITE(msgBuf,'(A30,1PE22.14)')
474     & ' TLM ref_cost_function =', fcref
475     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
476     & SQUEEZE_RIGHT, myThid )
477     WRITE(msgBuf,'(A30,1PE22.14)')
478     & ' TLM tangent-lin_grad =', ftlxxmemo
479     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
480     & SQUEEZE_RIGHT, myThid )
481     WRITE(msgBuf,'(A30,1PE22.14)')
482     & ' TLM finite-diff_grad =', gfd
483     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
484     & SQUEEZE_RIGHT, myThid )
485 jmc 1.32 #else
486 jmc 1.33 WRITE(msgBuf,'(A30,1PE22.14)')
487     & ' ADM ref_cost_function =', fcref
488     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
489     & SQUEEZE_RIGHT, myThid )
490     WRITE(msgBuf,'(A30,1PE22.14)')
491     & ' ADM adjoint_gradient =', adxxmemo
492     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
493     & SQUEEZE_RIGHT, myThid )
494     WRITE(msgBuf,'(A30,1PE22.14)')
495     & ' ADM finite-diff_grad =', gfd
496     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
497     & SQUEEZE_RIGHT, myThid )
498 heimbach 1.7 #endif
499    
500     print *, 'ph-grd ierr ---------------------------'
501     print *, 'ph-grd ierr = ', ierr, ', icomp = ', icomp,
502     & ', ichknum = ', ichknum
503    
504     c-- else of if ( ichknum ...
505     else
506     ierr_grdchk = -1
507    
508 jmc 1.31 c-- end of if ( ichknum ...
509 heimbach 1.2 endif
510    
511 jmc 1.31 c-- end of do icomp = ...
512 heimbach 1.7 enddo
513 heimbach 1.2
514 gforget 1.38 if (myProcId .EQ. grdchkwhichproc .AND. .NOT.useSingleCpuIO) then
515 jmc 1.31 CALL WRITE_REC_XYZ_RL(
516 heimbach 1.7 & 'grd_findiff' , tmpplot1, 1, 0, myThid)
517 jmc 1.31 CALL WRITE_REC_XYZ_RL(
518 heimbach 1.7 & 'grd_ratio_ad' , tmpplot2, 1, 0, myThid)
519 jmc 1.31 CALL WRITE_REC_XYZ_RL(
520 heimbach 1.7 & 'grd_ratio_ftl' , tmpplot3, 1, 0, myThid)
521     endif
522 heimbach 1.2
523 heimbach 1.7 c-- Everyone has to wait for the component to be reset.
524 jmc 1.37 c _BARRIER
525 heimbach 1.2
526     c-- Print the results of the gradient check.
527     call grdchk_print( ichknum, ierr_grdchk, mythid )
528    
529 heimbach 1.11 #endif /* ALLOW_GRDCHK */
530 heimbach 1.2
531 jmc 1.31 return
532 heimbach 1.2 end

  ViewVC Help
Powered by ViewVC 1.1.22