Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:45:04 UTC

view on githubraw file Latest commit b7b14afd on 2013-12-19 23:38:20 UTC
8702af1f36 Patr*0001       subroutine template()
                0002       use OAD_cp
                0003       use OAD_tape
                0004       use OAD_rev
                0005 
                0006 !$TEMPLATE_PRAGMA_DECLARATIONS
                0007 
                0008       integer :: cp_loop_variable_1,cp_loop_variable_2,
                0009      +     cp_loop_variable_3,cp_loop_variable_4
                0010 
                0011       type(modeType) :: our_orig_mode
                0012 
                0013       integer iaddr
                0014       external iaddr
                0015 
b7b14afd78 Jean*0016       INTEGER numItersHelper, myThidHelper, nIterMinHelper
8702af1f36 Patr*0017 
                0018       Real*8  cg2d_b_p(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0019       Real*8  cg2d_x_p(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0020 
                0021 #ifdef OAD_DEBUG_JOINT
                0022       character*(80):: indentation='                                        
                0023      +                                         '
                0024       our_indent=our_indent+1
                0025 
                0026       write(standardmessageunit, '(A,A,A)', ADVANCE='NO') 
                0027      +'OAD:',indentation(1:our_indent), 'enter __SRNAME__:'
                0028       call oad_dump_revmod(); call oad_dump_tapestats()
                0029       write(standardmessageunit,*) 
                0030 #endif
                0031       if (our_rev_mode%plain) then
                0032 #ifdef OAD_DEBUG_JOINT
                0033          write(standardmessageunit,'(A,A,A)') 
                0034      +'OAD:',indentation(1:our_indent), 
                0035      +' __SRNAME__: entering plain'
                0036 #endif
                0037 c set up for plain execution
                0038          our_orig_mode=our_rev_mode
                0039          our_rev_mode%arg_store=.FALSE.
                0040          our_rev_mode%arg_restore=.FALSE.
                0041          our_rev_mode%plain=.TRUE.
                0042          our_rev_mode%tape=.FALSE.
                0043          our_rev_mode%adjoint=.FALSE.
                0044 #ifdef OAD_DEBUG_JOINT
                0045          write(standardmessageunit,'(A,A,A)') 
                0046      +'OAD:',indentation(1:our_indent), 
                0047      +' __SRNAME__: runninng plain / down plain'
                0048 #endif
                0049          cg2d_b_p=cg2d_b%v
                0050          cg2d_x_p=cg2d_x%v
                0051         call cg2d (
                0052      I                cg2d_b_p,
                0053      U                cg2d_x_p,
                0054      O                firstResidual,
                0055      O                minResidualSq,
                0056      O                lastResidual,
                0057      U                numIters,
                0058      O                nIterMin,
                0059      I                myThid) 
                0060         cg2d_x%v=cg2d_x_p        
                0061 c reset the mode
                0062          our_rev_mode=our_orig_mode
                0063       end if
                0064       if (our_rev_mode%tape) then
                0065 
                0066          numItersHelper=numiters
                0067          mythidHelper=mythid
b7b14afd78 Jean*0068          nIterMinHelper=nIterMin
8702af1f36 Patr*0069 #ifdef OAD_DEBUG_JOINT
                0070          write(standardmessageunit,'(A,A,A)') 
                0071      +'OAD:',indentation(1:our_indent), 
                0072      +' __SRNAME__: entering tape'
                0073 #endif
                0074 c set up for plain execution
                0075          our_orig_mode=our_rev_mode
                0076          our_rev_mode%arg_store=.FALSE.
                0077          our_rev_mode%arg_restore=.FALSE.
                0078          our_rev_mode%plain=.TRUE.
                0079          our_rev_mode%tape=.FALSE.
                0080          our_rev_mode%adjoint=.FALSE.
                0081 #ifdef OAD_DEBUG_JOINT
                0082          write(standardmessageunit,'(A,A,A)') 
                0083      +'OAD:',indentation(1:our_indent), 
                0084      +' __SRNAME__: runninng plain / down plain'
                0085 #endif
                0086          cg2d_b_p=cg2d_b%v
                0087          cg2d_x_p=cg2d_x%v
                0088         call cg2d (
                0089      I                cg2d_b_p,
                0090      U                cg2d_x_p,
                0091      O                firstResidual,
                0092      O                minResidualSq,
                0093      O                lastResidual,
                0094      U                numIters,
                0095      O                nIterMin,
                0096      I                myThid) 
                0097         cg2d_x%v=cg2d_x_p        
                0098 c reset the mode
                0099          our_rev_mode=our_orig_mode
b7b14afd78 Jean*0100 c manually push integers to the tape:
                0101          call oad_tape_push(numItersHelper)
                0102          call oad_tape_push(nIterMinHelper)
                0103          call oad_tape_push(mythidHelper)
8702af1f36 Patr*0104       end if
                0105       if (our_rev_mode%adjoint) then
                0106 #ifdef OAD_DEBUG_JOINT
                0107          write(standardmessageunit,'(A,A,A)') 
                0108      +'OAD:',indentation(1:our_indent), 
                0109      +' __SRNAME__: entering adjoint'
                0110 #endif
                0111 c manually pop two integers from the tape:        
b7b14afd78 Jean*0112          call oad_tape_pop(mythid)
                0113          call oad_tape_pop(nIterMin)
                0114          call oad_tape_pop(numiters)
8702af1f36 Patr*0115 c selfadjoint:
                0116 c the original is called with 
                0117 c  cg2d(b,x,...)
                0118 c in the adjoint context if we 
                0119 c use the same code base 
                0120 c we call with 
                0121 c  cg2d(x_bar,bh,...
                0122 c where afterwards
                0123 c b_bar+=bh  and x_bar=0
                0124 c the adjoint second formal argument cg2d_x should be
                0125 c the values of the first argument:  
                0126       cg2d_b_p=cg2d_x%d
                0127 c the first formal argument cg2d_b should hold 
                0128 c the increment, i.e. we nullify the second formal 
                0129 c argument (cg2d_x) value:
                0130       cg2d_x_p=0.0
                0131 c set up for plain execution
                0132          our_orig_mode=our_rev_mode
                0133          our_rev_mode%arg_store=.FALSE.
                0134          our_rev_mode%arg_restore=.FALSE.
                0135          our_rev_mode%plain=.TRUE.
                0136          our_rev_mode%tape=.FALSE.
                0137          our_rev_mode%adjoint=.FALSE.
                0138 #ifdef OAD_DEBUG_JOINT
                0139          write(standardmessageunit,'(A,A,A)') 
                0140      +'OAD:',indentation(1:our_indent), 
                0141      +' __SRNAME__: runninng self adjoint / down plain'
                0142 #endif
                0143         call cg2d (
                0144      I                cg2d_b_p,
                0145      U                cg2d_x_p,
                0146      O                firstResidual,
                0147      O                minResidualSq,
                0148      O                lastResidual,
                0149      U                numIters,
                0150      O                nIterMin,
                0151      I                myThid) 
                0152 c reset the mode
                0153          our_rev_mode=our_orig_mode
                0154 c now the adjoint result is the increment 
                0155 c contained in the second formal argument
                0156          cg2d_b%d= cg2d_b%d+cg2d_x_p
                0157          cg2d_x%d=0.0
                0158       end if 
                0159 #ifdef OAD_DEBUG_JOINT
                0160       write(standardmessageunit,'(A,A,A)', ADVANCE='NO') 
                0161      +'OAD:',indentation(1:our_indent), 'leave __SRNAME__:'
                0162       call oad_dump_revmod(); call oad_dump_tapestats()
                0163       write(standardmessageunit,*) 
                0164 
                0165       our_indent=our_indent-1
                0166 #endif
                0167 
                0168       end subroutine template