Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:38:16 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
6d54cf9ca1 Ed H*0001 #include "CAL_OPTIONS.h"
a63ed37559 Patr*0002 
                0003       subroutine cal_StepsForDay(
                0004      I                            iday,
                0005      O                            firststep,
                0006      O                            laststep,
                0007      O                            nsteps,
                0008      I                            mythid
                0009      &                          )
                0010 
                0011 c     ==================================================================
                0012 c     SUBROUTINE cal_StepsForDay
                0013 c     ==================================================================
                0014 c
                0015 c     o Given the current day of the integration this routine returns
                0016 c       first, the last and the number of model steps the will have to
                0017 c       be performed.
                0018 c
                0019 c       This routine also checks consistency of variables quite
                0020 c       extensively.
                0021 c
                0022 c     started: Christian Eckert eckert@mit.edu  06-Apr-2000
                0023 c
c3cd6c250f Jean*0024 c     changed:
a63ed37559 Patr*0025 c
                0026 c     ==================================================================
                0027 c     SUBROUTINE cal_StepsForDay
                0028 c     ==================================================================
                0029 
                0030       implicit none
                0031 
                0032 c     == global variables ==
                0033 
                0034 #include "cal.h"
                0035 
                0036 c     == routine arguments ==
                0037 
                0038       integer iday
                0039       integer firststep
                0040       integer laststep
                0041       integer nsteps
                0042       integer mythid
                0043 
                0044 c     == local variables ==
                0045 
                0046       integer ierr
                0047       integer mdstep
                0048       integer numdays
                0049       integer numsteps
                0050       integer frac1
                0051       integer frac2
                0052       integer frac3
                0053       integer frac4
                0054       integer fullsteps
                0055       integer firstyear
                0056       integer firstmonth
                0057       integer firstday
                0058       integer lyfirst
                0059       integer startsecs
                0060       integer lastday
                0061       integer endsecs
                0062 
                0063 c     == external ==
                0064 
                0065       external cal_IntDays
                0066       integer  cal_IntDays
                0067       external cal_IsLeap
                0068       integer  cal_IsLeap
                0069 
                0070 c     == end of interface ==
                0071 
                0072       numdays    = cal_IntDays( mythid )
                0073       lyfirst    = cal_IsLeap( firstyear, mythid )
                0074 
                0075       mdstep     = int(modelstep)
                0076 
                0077       firstyear  = modelstartdate(1)/10000
                0078       firstmonth = mod(modelstartdate(1)/100,100)
                0079       firstday   = mod(modelstartdate(1),100)
                0080       lastday    = mod(modelenddate(1),100)
                0081 
                0082       startsecs  = (modelstartdate(2)/10000)*secondsperhour +
                0083      &             mod(modelstartdate(2)/100,100)*secondsperminute +
                0084      &             mod(modelstartdate(2),100)
                0085       endsecs    = (modelenddate(2)/10000)*secondsperhour +
                0086      &             mod(modelenddate(2)/100,100)*secondsperminute +
                0087      &             mod(modelenddate(2),100)
                0088 
                0089       if ( numdays .eq. 1 ) then
                0090         if ( iday .eq. firstday ) then
                0091 c--       Get the number of steps in the first day.
                0092           if ( firstday .eq. lastday ) then
                0093             firststep = 1
                0094             laststep  = modelintsteps
                0095           else if ( mod(firstday+1,ndaymonth(firstmonth,lyfirst)) .eq.
                0096      &              lastday ) then
                0097 c--         This can only happen if we end at midnight of the next day.
                0098             if ( modelenddate(2) .eq. 0 ) then
                0099               firststep = 1
                0100               laststep  = modelintsteps
                0101 c--           Note: This holds only if modelenddate was determined
c3cd6c250f Jean*0102 c--                 such that it coincides with the model final time.
a63ed37559 Patr*0103             else
                0104 c--           We do not end at midnight of the first day of
                0105 c--           the next month.
                0106               ierr = 2604
                0107               call cal_PrintError( ierr, mythid )
                0108               stop ' stopped in cal_StepsForDay.'
                0109             endif
                0110           else
                0111 c--         The first and the last day are inconsistent with iday.
                0112             ierr = 2603
                0113             call cal_PrintError( ierr, mythid )
                0114             stop ' stopped in cal_StepsForDay.'
                0115           endif
                0116         else
                0117 c--       The variables numdays and iday are inconsistent;
                0118 c--       ( iday .gt. numdays ).
                0119           ierr = 2602
                0120           call cal_PrintError( ierr, mythid )
                0121           stop ' stopped in cal_StepsForDay.'
                0122         endif
                0123 
                0124       else if ( numdays .gt. 1 ) then
                0125 c--     More than one day of integration.
                0126         if ( iday .eq. 1 ) then
                0127           firststep = 1
                0128           laststep  = int((secondsperday - startsecs)/mdstep)
                0129         else if ( ( iday .gt. 1      ) .and.
                0130      &            ( iday .lt. numdays) ) then
                0131 c--       Somewhere between first and last month.
                0132 c--       The first steps in iday.
                0133           fullsteps = int((secondsperday - startsecs)/mdstep)
                0134           numsteps  = fullsteps
c3cd6c250f Jean*0135 c--       What is left in the first day (frac1).
a63ed37559 Patr*0136           frac1     = (secondsperday - startsecs) - fullsteps*mdstep
                0137           fullsteps = int(secondsperday/modelstep)
c3cd6c250f Jean*0138 c--       What is left in a complete day (frac2).
a63ed37559 Patr*0139           frac2     = secondsperday - fullsteps*mdstep
c3cd6c250f Jean*0140 c--       What is left up to the current day (frac3).
a63ed37559 Patr*0141           frac3     = frac1 + frac2*(iday - 1)
                0142           numsteps  = numsteps + (iday - 1)*fullsteps +
                0143      &                frac3/mdstep
                0144           laststep  = numsteps
                0145           firststep = laststep - secondsperday/mdstep + 1
                0146 
                0147         else if ( iday .eq. numdays ) then
                0148 c--       The last day of integration.
                0149 c--       The first step in iday.
                0150           fullsteps = int((secondsperday - startsecs)/mdstep)
                0151           numsteps  = fullsteps
c3cd6c250f Jean*0152 c--       What is left in the first day (frac1).
a63ed37559 Patr*0153           frac1     = (secondsperday - startsecs) - fullsteps*mdstep
                0154           fullsteps = int(secondsperday/modelstep)
c3cd6c250f Jean*0155 c--       What is left in a complete day (frac2).
a63ed37559 Patr*0156           frac2     = secondsperday - fullsteps*mdstep
c3cd6c250f Jean*0157 c--       What is left up to the day before the last (frac3).
a63ed37559 Patr*0158           frac3     = frac1 + frac2*(iday - 2)
                0159           numsteps  = numsteps + (iday - 2)*fullsteps
                0160 c--       The last step in iday.
                0161           if ( modelenddate(2) .eq. 0 ) then
                0162 c--         This can only happen if we end at midnight of the next day.
                0163             laststep  = numsteps + fullsteps
                0164             firststep = numsteps + 1
                0165 c--         Note: There should be no fraction left
                0166 c--               ( mod(frac3,mdstep) = frac3/mdstep ) if modelenddate
                0167 c--               is based on an integral number of timesteps.
                0168           else
                0169             frac4     = frac3 + endsecs
                0170             numsteps  = numsteps + frac4/mdstep
                0171             laststep  = numsteps
                0172 
                0173 c--         Note: There should be no fraction left
                0174 c--               ( mod(frac4,mdstep = frac4/mdstep ) if modelenddate
                0175 c--               is based on an integral number of timesteps.
                0176 
                0177             firststep = laststep - endsecs/mdstep + 1
                0178           endif
                0179         else
                0180 c--       The variables iday and numdays are inconsistent.
                0181           ierr = 2605
                0182           call cal_PrintError( ierr, mythid )
                0183           stop ' stopped in cal_DaysForMonth.'
                0184         endif
                0185       else
                0186 c--     The number of days to integrate is wrong; check cal_IntDays.
                0187         ierr = 2601
                0188         call cal_PrintError( ierr, mythid )
                0189         stop ' stopped in cal_StepsForDay.'
                0190       endif
                0191 
                0192 c--   The number of days to integrate in the given month.
                0193       nsteps = laststep - firststep + 1
                0194 
                0195       return
                0196       end
                0197