/[MITgcm]/MITgcm/tools/OAD_support/ad_template.revolve.F
ViewVC logotype

Annotation of /MITgcm/tools/OAD_support/ad_template.revolve.F

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


Revision 1.4 - (hide annotations) (download)
Fri Jun 21 19:58:16 2013 UTC (10 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint64o, checkpoint64q, checkpoint64p, checkpoint64r, checkpoint64n, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l
Changes since 1.3: +6 -0 lines
OAD streamice support

1 heimbach 1.3 #include "PACKAGES_CONFIG.h"
2    
3 utke 1.1 subroutine template()
4     use OAD_cp
5     use OAD_tape
6     use OAD_rev
7     use revolve
8    
9     c we may need these for the checkpointing
10     use SIZE_mod
11     use EEPARAMS_mod
12     use PARAMS_mod
13     use BAR2_mod
14     use BARRIER_mod
15 heimbach 1.3 #ifdef ALLOW_CD_CODE
16 utke 1.1 use CD_CODE_VARS_mod
17 heimbach 1.3 #endif
18 utke 1.1 use CG2D_mod
19     use CG3D_mod
20     use DYNVARS_mod
21     use EESUPPORT_mod
22     use EOS_mod
23     use EXCH_mod
24     use FC_NAMEMANGLE_mod
25     use FFIELDS_mod
26     use GAD_mod
27     use GLOBAL_MAX_mod
28     use GLOBAL_SUM_mod
29 heimbach 1.3 #ifdef ALLOW_GMREDI
30 utke 1.1 use GMREDI_mod
31     use GMREDI_TAVE_mod
32 heimbach 1.3 #endif
33 utke 1.1 use GRID_mod
34 heimbach 1.3 use MOM_VISC_mod
35 utke 1.1 use MPI_INFO_mod
36 heimbach 1.3 #ifdef ALLOW_SHAP_FILT
37     use SHAP_FILT_mod
38     #endif
39 heimbach 1.4 #ifdef ALLOW_STREAMICE
40     use STREAMICE_mod
41     use STREAMICE_ADV_mod
42     use STREAMICE_BDRY_mod
43     use STREAMICE_CG_mod
44     #endif
45 utke 1.1 use SOLVE_FOR_PRESSURE3D_mod
46     use SOLVE_FOR_PRESSURE_mod
47     use SURFACE_mod
48     use tamc_mod
49     use tamc_keys_mod
50     use cost_mod
51     use g_cost_mod
52     use ctrl_mod
53     use ctrl_dummy_mod
54     use ctrl_weights_mod
55     use optim_mod
56     use grdchk_mod
57    
58     !$TEMPLATE_PRAGMA_DECLARATIONS
59     LOGICAL :: initialized=.FALSE.
60     TYPE(rvAction),save :: theAction
61 utke 1.2 CHARACTER(80) :: errorMsg
62 utke 1.1 integer, save :: jointCPCount
63     integer, save :: currIter
64    
65     integer :: cp_loop_variable_1,cp_loop_variable_2,
66     + cp_loop_variable_3,cp_loop_variable_4,cp_loop_variable_5
67    
68     type(modeType) :: our_orig_mode
69    
70     integer iaddr
71     external iaddr
72    
73     #ifdef OAD_DEBUG_JOINT
74     character*(80):: indentation='
75     + '
76     our_indent=our_indent+1
77    
78     write(standardmessageunit, '(A,A,A)', ADVANCE='NO')
79     +'OAD:',indentation(1:our_indent), 'enter __SRNAME__:'
80     call oad_dump_revmod(); call oad_dump_tapestats()
81     write(standardmessageunit,*)
82     #endif
83    
84     nIter0 = NINT( (startTime-baseTime)/deltaTClock )
85     if (our_rev_mode%arg_store) then
86     call cp_write_open()
87     #ifdef OAD_DEBUG_JOINT
88     write(standardmessageunit,'(A,A,A)')
89     +'OAD:',indentation(1:our_indent),
90     +' __SRNAME__: entering arg store'
91     #endif
92     !$PLACEHOLDER_PRAGMA$ id=8
93     call cp_close()
94     end if
95     if (our_rev_mode%arg_restore) then
96     #ifdef OAD_DEBUG_JOINT
97     write(standardmessageunit,'(A,A,A)')
98     +'OAD:',indentation(1:our_indent),
99     +' __SRNAME__: entering arg restore'
100     #endif
101     call cp_read_open()
102     !$PLACEHOLDER_PRAGMA$ id=9
103     call cp_close()
104     end if
105     if (our_rev_mode%plain) then
106     #ifdef OAD_DEBUG_JOINT
107     write(standardmessageunit,'(A,A,A)')
108     +'OAD:',indentation(1:our_indent),
109     +' __SRNAME__: run plain, down plain'
110     #endif
111     DO iloop = 1, nTimeSteps
112     CALL OpenAD_forward_step( iloop, mytime, myiter, mythid )
113     enddo
114     end if
115     if (our_rev_mode%tape) then
116     #ifdef OAD_DEBUG_JOINT
117     write(standardmessageunit,'(A,A,A)')
118     +'OAD:',indentation(1:our_indent),
119     +' __SRNAME__: run tape, down revolve until first U turn'
120     #endif
121     currIter=0
122     jointCPcount=cp_fNumber()
123     initialized=rvInit(nTimeSteps,120,
124     + errorMsg,theAction)
125     IF (.NOT.initialized) WRITE(*,'(A,A)') 'Error: ', errorMsg
126     do while (theAction%actionFlag/=rvDone)
127     theAction=rvNextAction()
128     select case (theAction%actionFlag)
129     case (rvStore)
130     call cp_write_open(theAction%cpNum+jointCPCount)
131     !$PLACEHOLDER_PRAGMA$ id=8
132     call cp_close
133     case (rvForward)
134     call OAD_revPlain
135     do currIter=currIter,theAction%iteration-1
136     CALL OpenAD_forward_step( currIter+1, mytime,
137     +myiter, mythid )
138     end do
139     call OAD_revTape
140     case (rvFirstUTurn)
141     CALL OpenAD_forward_step( currIter+1, mytime, myiter,
142     +mythid )
143     ! get out now ...
144     exit
145     end select
146     end do
147     end if
148     if (our_rev_mode%adjoint) then
149     IF (.NOT.initialized) WRITE(*,'(A)') 'Error: not initialized'
150     do while (theAction%actionFlag/=rvDone)
151     select case (theAction%actionFlag)
152     case (rvFirstUTurn)
153     !we taped already ... see above
154     CALL OpenAD_forward_step( currIter+1, mytime, myiter,
155     +mythid )
156     case (rvStore)
157     call cp_write_open(theAction%cpNum+jointCPCount)
158     !$PLACEHOLDER_PRAGMA$ id=8
159     call cp_close
160     case (rvRestore)
161     call cp_read_open(theAction%cpNum+jointCPCount)
162     !$PLACEHOLDER_PRAGMA$ id=9
163     currIter=theAction%iteration
164     call cp_close
165     case (rvForward)
166     call OAD_revPlain
167     do currIter=currIter,theAction%iteration-1
168     CALL OpenAD_forward_step( currIter+1, mytime, myiter,
169     + mythid )
170     end do
171     call OAD_revAdjoint
172     case (rvUTurn)
173     call OAD_revTape
174     CALL OpenAD_forward_step( currIter+1, mytime, myiter,
175     +mythid )
176     call OAD_revAdjoint
177     CALL OpenAD_forward_step( currIter+1, mytime, myiter,
178     +mythid )
179     end select
180     theAction=rvNextAction()
181     end do
182     end if
183    
184     #ifdef OAD_DEBUG_JOINT
185     write(standardmessageunit,'(A,A,A)', ADVANCE='NO')
186     +'OAD:',indentation(1:our_indent), 'leave __SRNAME__:'
187     call oad_dump_revmod(); call oad_dump_tapestats()
188     write(standardmessageunit,*)
189    
190     our_indent=our_indent-1
191     #endif
192    
193     end subroutine template

  ViewVC Help
Powered by ViewVC 1.1.22