/[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.32 - (hide annotations) (download)
Fri Aug 12 18:36:49 2011 UTC (12 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63b
Changes since 1.31: +26 -30 lines
send adjoint gradient to all procs so that all procs print it to STDOUT
 (a hack to get gradient-check tested with testreport)

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

  ViewVC Help
Powered by ViewVC 1.1.22