psmile_check_action.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_Create_timeaxis
00008 !
00009 ! !INTERFACE:
00010 !
00011 subroutine psmile_check_action ( field_id, task_id, request, &
00012                                  julian_day, julian_dayb,    &
00013                                  julian_sec, julian_secb,    &
00014                                  action )
00015 !
00016 ! !USES:
00017 !
00018   use PSMILe, dummy_interface => psmile_check_action
00019   use PSMILe_SMIOC, only : sga_smioc_comp, transient
00020   use PRISM_constants
00021 
00022   implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026   Integer, Intent (In)  :: field_id
00027   Integer, Intent (In)  :: task_id
00028 
00029   Double Precision, Intent (In) :: julian_day, julian_dayb(2)
00030   Double Precision, Intent (In) :: julian_sec, julian_secb(2)
00031 
00032   Logical, Intent (In)  :: request    ! Switch to distinguish between calls
00033                                       ! from prism_put_inquire and prism_put
00034 !
00035 ! !OUTPUT PARAMETERS:
00036 !
00037   Logical, Intent (Out) :: action(3)  ! 1: coupling action
00038                                       ! 2: IO action
00039                                       ! 3: restart action
00040 !
00041 ! !LOCAL VARIABLES
00042 !
00043   Double Precision      :: upper_bnd, lower_bnd, test_bnd
00044   Double Precision      :: juday_axis
00045 
00046   Integer               :: i, lag
00047   Logical               :: time_action 
00048 
00049   Type (transient), pointer :: sga_smioc_transi(:)
00050 !
00051 ! !DESCRIPTION:
00052 !
00053 !    Subroutine "psmile_check_action" determines whether a prism_put or
00054 !    prism_get call has to activated based on the data information and
00055 !    the frequency of coupling or IO events as it is requested in the
00056 !    SMIOC file.
00057 !
00058 ! !REVISION HISTORY:
00059 !   Date      Programmer   Description
00060 ! ----------  ----------   -----------
00061 ! 04.01.26    R. Redler    created
00062 ! 11.02.24    R. Redler    revised
00063 !
00064 !EOP
00065 !------------------------------------------------------------------------
00066 !
00067 ! $Id: psmile_check_action.F90 3248 2011-06-23 13:03:19Z coquart $
00068 ! $Author: coquart $
00069 !
00070   Character(len=len_cvs_string), save :: mycvs = 
00071       '$Id: psmile_check_action.F90 3248 2011-06-23 13:03:19Z coquart $'
00072 !
00073 !------------------------------------------------------------------------
00074 
00075 #ifdef VERBOSE
00076   print *, trim(ch_id), ': psmile_check_action: field_id ', field_id
00077 
00078   call psmile_flushstd
00079 #endif /* VERBOSE */
00080   
00081   !----------------------------------------------------------------------
00082   ! 1st   Initialisation
00083   !----------------------------------------------------------------------
00084   
00085   sga_smioc_transi => sga_smioc_comp(Fields(field_id)%comp_id)%sga_smioc_transi
00086 
00087   time_action = .false.
00088 
00089   lower_bnd  = julian_dayb(1) + julian_secb(1)/86400.0
00090   upper_bnd  = julian_dayb(2) + julian_secb(2)/86400.0
00091 
00092 #ifdef DEBUG
00093   print *, trim(ch_id), ': psmile_check_action: julian_dayb(1), julian_secb(1)', &
00094                                                 julian_dayb(1), julian_secb(1)
00095   print *, trim(ch_id), ': psmile_check_action: julian_dayb(2), julian_secb(2)', &
00096                                                 julian_dayb(2), julian_secb(2)
00097   call psmile_flushstd
00098 #endif
00099 
00100   action = .false.
00101 
00102   !----------------------------------------------------------------------
00103   ! 2nd   Input section (PRISM_get)
00104   !----------------------------------------------------------------------
00105 
00106   if ( task_id == 0 ) then
00107 
00108      !-------------------------------------------------------------------
00109      ! ...   Check date and time
00110      !-------------------------------------------------------------------
00111 
00112 #ifdef DEBUG
00113      print *, trim(ch_id), ': psmile_check_action: Event at ',    &
00114                        Fields(field_id)%Taskin%Judate_Event%days, &
00115                        Fields(field_id)%Taskin%Judate_Event%secs
00116      call psmile_flushstd
00117 #endif
00118 
00119      if ( lower_bnd > Fields(field_id)%Taskin%Judate_Stop%days &
00120                     + Fields(field_id)%Taskin%Judate_Stop%secs/86400.0 ) then
00121         time_action = .false.
00122      endif
00123 
00124      if ( upper_bnd < Fields(field_id)%Taskin%Judate_Start%days &
00125                     + Fields(field_id)%Taskin%Judate_Start%secs/86400.0 ) then
00126         time_action = .false.
00127      endif
00128 
00129      !----------------------------------------------------------------
00130      ! ...   Update Event if neccessary
00131      !----------------------------------------------------------------
00132 
00133      test_bnd = Fields(field_id)%Taskin%Judate_Event%days &
00134               + Fields(field_id)%Taskin%Judate_Event%secs/86400.0
00135 
00136      do while ( lower_bnd >= test_bnd )
00137 
00138         call psmile_calc_new_date ( Fields(field_id)%Taskin%Judate_Event, &
00139                                     Fields(field_id)%Taskin%delta_time )
00140 #ifdef DEBUG
00141         print *, trim(ch_id), ': psmile_check_action: Event shifted to ', &
00142              Fields(field_id)%Taskin%Judate_Event%days, &
00143              Fields(field_id)%Taskin%Judate_Event%secs
00144         call psmile_flushstd
00145 #endif
00146         test_bnd = Fields(field_id)%Taskin%Judate_Event%days &
00147                  + Fields(field_id)%Taskin%Judate_Event%secs/86400.0
00148      enddo
00149 
00150      juday_axis = Fields(field_id)%Taskin%Judate_Event%days &
00151                 + Fields(field_id)%Taskin%Judate_Event%secs/86400.0
00152 
00153      if ( lower_bnd < juday_axis .and. upper_bnd >= juday_axis ) then
00154 
00155         time_action = .true.
00156 
00157         !-------------------------------------------------------------
00158         ! ...   Set next event if neccessary
00159         !-------------------------------------------------------------
00160 
00161         call psmile_calc_new_date ( Fields(field_id)%Taskin%Judate_Event, &
00162                                     Fields(field_id)%Taskin%delta_time ) 
00163 #ifdef DEBUG
00164         print *, trim(ch_id), ': psmile_check_action: Event shifted after action to ', &
00165                                             Fields(field_id)%Taskin%Judate_Event%days, &
00166                                             Fields(field_id)%Taskin%Judate_Event%secs
00167         call psmile_flushstd
00168 #endif
00169 
00170      endif
00171 
00172      !----------------------------------------------------------------
00173      ! ...   Determine input path
00174      !----------------------------------------------------------------
00175 
00176      if ( time_action ) then
00177 
00178         do i = 1, sga_smioc_transi(Fields(field_id)%smioc_loc)%sg_transi_in%ig_nb_in_orig
00179 
00180            if ( sga_smioc_transi(Fields(field_id)%smioc_loc)%sg_transi_in%sga_in_orig(i)%ig_orig_type == PSMILe_comp ) &
00181                 action(1) = .true.
00182 
00183            if ( sga_smioc_transi(Fields(field_id)%smioc_loc)%sg_transi_in%sga_in_orig(i)%ig_orig_type == PSMILe_file ) &
00184                 action(2) = .true.
00185 
00186         enddo
00187 
00188      endif
00189 
00190   endif ! task_id
00191 
00192   !----------------------------------------------------------------------
00193   ! 3rd   Output section (PRISM_put)
00194   !----------------------------------------------------------------------
00195 
00196   if ( task_id > 0 ) then
00197 
00198 #ifdef DEBUG
00199      print *, trim(ch_id), ': psmile_check_action: Event at ',              &
00200           Fields(field_id)%Taskout(task_id)%Judate_Event%days, &
00201           Fields(field_id)%Taskout(task_id)%Judate_Event%secs
00202      call psmile_flushstd
00203 #endif
00204 
00205      lag = sga_smioc_transi(Fields(field_id)%smioc_loc)%sga_transi_out(task_id)%ig_lag
00206 
00207      !-------------------------------------------------------------------
00208      ! ...   Check date and time
00209      !-------------------------------------------------------------------
00210 
00211      if ( lag == 0 .and. lag /= PSMILe_undef ) then
00212 
00213         if ( upper_bnd >= Fields(field_id)%Taskout(task_id)%Judate_Stop%days           &
00214                         + Fields(field_id)%Taskout(task_id)%Judate_Stop%secs/86400.0 ) &
00215              action(3) = .true.
00216 
00217      else if ( lag > 0 .and. lag /= PSMILe_undef ) then
00218 
00219         if ( upper_bnd > Fields(field_id)%Taskout(task_id)%Judate_Stop%days           &
00220                        + Fields(field_id)%Taskout(task_id)%Judate_Stop%secs/86400.0 ) &
00221              action(3) = .true.
00222      endif
00223 
00224      if ( lower_bnd > Fields(field_id)%Taskout(task_id)%Judate_Stop%days           &
00225                     + Fields(field_id)%Taskout(task_id)%Judate_Stop%secs/86400.0 ) &
00226           time_action = .false.
00227 
00228      if ( upper_bnd < Fields(field_id)%Taskout(task_id)%Judate_Start%days           &
00229                     + Fields(field_id)%Taskout(task_id)%Judate_Start%secs/86400.0 ) &
00230           time_action = .false.
00231 
00232      !----------------------------------------------------------------
00233      ! ...   Update Event if neccessary
00234      !----------------------------------------------------------------
00235 
00236      if ( request .or. .not. Fields(field_id)%Taskout(task_id)%requested ) then
00237 
00238         test_bnd = Fields(field_id)%Taskout(task_id)%Judate_Event%days &
00239                  + Fields(field_id)%Taskout(task_id)%Judate_Event%secs/86400.0
00240 
00241         do while ( lower_bnd >= test_bnd )
00242 
00243            call psmile_calc_new_date ( Fields(field_id)%Taskout(task_id)%Judate_Event, &
00244                                        Fields(field_id)%Taskout(task_id)%delta_time ) 
00245 #ifdef DEBUG
00246            print *, trim(ch_id), ': psmile_check_action: Event shifted to ', &
00247                         Fields(field_id)%Taskout(task_id)%Judate_Event%days, &
00248                         Fields(field_id)%Taskout(task_id)%Judate_Event%secs
00249            call psmile_flushstd
00250 #endif
00251            test_bnd = Fields(field_id)%Taskout(task_id)%Judate_Event%days &
00252                     + Fields(field_id)%Taskout(task_id)%Judate_Event%secs/86400.0
00253         enddo
00254 
00255      endif
00256 
00257      juday_axis = Fields(field_id)%Taskout(task_id)%Judate_Event%days &
00258                 + Fields(field_id)%Taskout(task_id)%Judate_Event%secs/86400.0
00259 
00260      !----------------------------------------------------------------
00261      ! ...   Check for an event
00262      !----------------------------------------------------------------
00263 
00264      if ( lower_bnd < juday_axis .and. upper_bnd >= juday_axis ) then
00265 
00266         time_action = .true.
00267 
00268         !-------------------------------------------------------------
00269         ! ...   Set next event if neccessary but only when called from
00270         !       an active prism_put 
00271         !-------------------------------------------------------------
00272 
00273         if ( .not. request ) then
00274 
00275            call psmile_calc_new_date ( Fields(field_id)%Taskout(task_id)%Judate_Event, &
00276                                        Fields(field_id)%Taskout(task_id)%delta_time ) 
00277 #ifdef DEBUG
00278            print *, trim(ch_id), ': psmile_check_action: Event shifted after action to ', &
00279                                      Fields(field_id)%Taskout(task_id)%Judate_Event%days, &
00280                                      Fields(field_id)%Taskout(task_id)%Judate_Event%secs
00281            call psmile_flushstd
00282 #endif
00283         endif
00284 
00285      endif
00286 
00287      !----------------------------------------------------------------
00288      ! ...   Determine output path
00289      !----------------------------------------------------------------
00290 
00291      if ( time_action ) then
00292 
00293         if ( sga_smioc_transi(Fields(field_id)%smioc_loc)%sga_transi_out(task_id)%ig_dest_type == PSMILe_comp ) &
00294              action(1) = .true.
00295 
00296         if ( sga_smioc_transi(Fields(field_id)%smioc_loc)%sga_transi_out(task_id)%ig_dest_type == PSMILe_file ) &
00297              action(2) = .true.
00298      endif
00299 
00300      ! For lags > 0 we can either write a couling restart for a particular
00301      ! task or send it to a transformer or application process, but not 
00302      ! both at the same time. For lag == 0 suppress writing of restart file
00303      ! at coupling steps. 
00304 
00305      if ( lag == 0 .and. action(1) ) action(3) = .false.
00306      if ( lag  > 0 .and. action(3) ) action(1) = .false.
00307 
00308      Fields(field_id)%Taskout(task_id)%requested = request
00309 
00310   endif ! task_id
00311 
00312 #ifdef VERBOSE
00313   print *, trim(ch_id), ': psmile_check_action: ', 'Cpl=',  action(1), &
00314        'IO=',   action(2), &
00315        'Rest=', action(3)
00316   print *, trim(ch_id), ': psmile_check_action: end'
00317   call psmile_flushstd
00318 #endif /* VERBOSE */
00319 
00320 end subroutine psmile_check_action
00321 
00322 ! =======================================================================
00323 
00324 subroutine psmile_calc_new_date ( date, delta_time )
00325   !
00326   ! !USES:
00327   !
00328   use PRISM_constants, ONLY : PRISM_Time_Struct
00329   use PSMILe, ONLY : PSMILe_Time_Struct
00330   use PRISM_calendar, ONLY : psmile_date2ju, psmile_ju2date
00331   !
00332   ! !INPUT/OUTOUT PARAMETERS:
00333   !
00334   Type(PSMILe_Time_Struct), INTENT(INOUT) :: date
00335   !
00336   ! !INPUT PARAMETERS:
00337   !
00338   Type(PRISM_Time_Struct), INTENT(IN)     :: delta_time
00339   !
00340   ! !LOCAL VARIABLES
00341   !
00342   Type(PRISM_Time_Struct ) :: tmp_date(2)
00343 
00344   Double Precision   :: temp_days
00345   Double Precision   :: temp_secs, secs
00346 
00347   Integer, Parameter :: iold = 1
00348   Integer, Parameter :: inew = 2
00349 
00350   Integer            :: add_days
00351   !
00352   ! !DESCRIPTION:
00353   !
00354   !    Subroutine "psmile_calc_new_date" updates the date/time
00355   !      provided as Julian data/time with the delta time provided
00356   !      by the user in the XML file.
00357   !
00358   ! !REVISION HISTORY:
00359   !   Date      Programmer   Description
00360   ! ----------  ----------   -----------
00361   ! 11.02.24    R. Redler    created
00362   !
00363   !EOP
00364   !----------------------------------------------------------------------
00365 
00366   !----------------------------------------------------------------------
00367   ! Convert Julian date/time into calendar
00368   !----------------------------------------------------------------------
00369 
00370   call psmile_ju2date ( tmp_date(iold), date%days, date%secs )
00371 
00372   !----------------------------------------------------------------------
00373   ! Update of years and months first
00374   !----------------------------------------------------------------------
00375 
00376   tmp_date(inew)%year  = tmp_date(iold)%year  + delta_time%year
00377   tmp_date(inew)%month = tmp_date(iold)%month + delta_time%month
00378 
00379   if ( tmp_date(inew)%month > 12 ) then
00380        tmp_date(inew)%month = tmp_date(inew)%month - 12
00381        tmp_date(inew)%year  = tmp_date(inew)%year  +  1
00382   endif
00383 
00384   !----------------------------------------------------------------------
00385   ! Updates of time in Julian time
00386   !----------------------------------------------------------------------
00387 
00388   tmp_date(inew)%day    = tmp_date(iold)%day
00389   tmp_date(inew)%hour   = tmp_date(iold)%hour
00390   tmp_date(inew)%minute = tmp_date(iold)%minute
00391   tmp_date(inew)%second = tmp_date(iold)%second
00392 
00393   call psmile_date2ju ( tmp_date(inew), temp_days, temp_secs )
00394 
00395   secs = delta_time%hour*3600.0 + delta_time%minute*60.0 &
00396        + delta_time%second      + temp_secs
00397 
00398   add_days = int(secs/86400.0)
00399 
00400   !----------------------------------------------------------------------
00401   ! Update input Julian date by adding days and seconds
00402   !----------------------------------------------------------------------
00403 
00404   date%days = temp_days + delta_time%day + add_days
00405   date%secs = secs - (float(add_days) * 86400.0)
00406 
00407 end subroutine psmile_calc_new_date

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1