psmile_get_restart.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PSMILe_Get_restart
00008 !
00009 ! !INTERFACE:
00010 
00011 subroutine psmile_get_restart ( ierror )
00012 !
00013 ! !USES:
00014 !
00015   use PRISM_constants
00016   use PRISM_calendar
00017 
00018   use PSMILe
00019   use PSMILe_SMIOC, only : sga_smioc_comp, transient
00020 
00021   implicit none
00022 !
00023 ! !OUTPUT PARAMETERS:
00024 !
00025   Integer, Intent (Out)           :: ierror
00026 !
00027 !     Returns the error code of psmile_get_restart;
00028 !             ierror = 0 : No error
00029 !             ierror > 0 : Severe error
00030 !
00031 !
00032 ! !LOCAL VARIABLES
00033 !
00034   Integer    :: comp_id
00035   Integer    :: field_id
00036   Integer    :: task_id
00037   Integer    :: smioc_loc     ! reference to location in SMIOC structure
00038   Integer    :: io_taskid     ! task Id shifted into range of restarts
00039   Integer    :: nb_transiouts ! number of transients for each task
00040   Integer    :: nbr_fields_of_comp ! number of fields defined per component
00041 
00042   Integer    :: nbr_fields    ! number of fields in a bundle
00043   Integer    :: datatype
00044   Integer    :: lag           ! time lag set in SMIOC
00045   Integer    :: len           ! length of data array
00046   Integer    :: info          ! 1 if valid restart was found, 0 otherwise
00047 
00048   Integer, Parameter :: nerrp = 2
00049   Integer            :: ierrp (nerrp)
00050 
00051   Logical            :: local_timeop
00052 
00053   Integer,          Allocatable :: idata_array(:)
00054   Real,             Allocatable :: rdata_array(:)
00055   Double Precision, Allocatable :: ddata_array(:)
00056 
00057   Double Precision   :: julian_day, julian_dayb(2)
00058   Double Precision   :: julian_sec, julian_secb(2)
00059 
00060   Type (GridFunction), Pointer :: fp
00061   Type (transient), Pointer    :: sga_smioc_transi(:)
00062 !
00063 ! !DESCRIPTION:
00064 !
00065 ! Subroutine "psmile_get_restart" is called by prism_enddef and reads
00066 !    the data from a restart file; necessary only when lag was set
00067 !    to something different from zero. Restart data are then sent to
00068 !    the transformer or directly to the applications via psmile_put_field_
00069 !
00070 ! !REVISION HISTORY:
00071 !
00072 !   Date      Programmer   Description
00073 ! ----------  ----------   -----------
00074 ! 04.08.09    R. Redler    created
00075 ! 04.08.29    R. Vogelsang psmile_check_restart is called with the unshifted
00076 !                          taskid.
00077 !
00078 !EOP
00079 !----------------------------------------------------------------------
00080 !
00081 ! $Id: psmile_get_restart.F90 2687 2010-10-28 15:15:52Z coquart $
00082 ! $Autor$
00083 !
00084   Character(len=len_cvs_string), save :: mycvs = 
00085        '$Id: psmile_get_restart.F90 2687 2010-10-28 15:15:52Z coquart $'
00086 !
00087 !----------------------------------------------------------------------
00088 !
00089 #ifdef VERBOSE
00090   print 9980, trim(ch_id)
00091   call psmile_flushstd
00092 #endif /* VERBOSE */
00093 
00094 !-----------------------------------------------------------------------
00095 ! 1st Initialization
00096 !-----------------------------------------------------------------------
00097 
00098   ierror = 0
00099   nb_transiouts = 0
00100   local_timeop = .false.
00101 
00102 !-----------------------------------------------------------------------
00103 ! ... set julian date and bounds for reading the restart. The restart
00104 !     date should be equal to the job start date. Nevertheless we allow
00105 !     for an uncertainty of 2 seconds. Following the ISO 8601 standard
00106 !     we should not allow times like 24:00 but used 00:00 instead.
00107 !-----------------------------------------------------------------------
00108 
00109   call psmile_date2ju(PRISM_Jobstart_date, julian_day, julian_sec)
00110 
00111   if ( julian_sec >= 86398.0 ) then
00112      julian_dayb(1) = julian_day
00113      julian_dayb(2) = julian_day + 1.0
00114      julian_secb(1) = julian_sec - 2.0
00115      julian_secb(2) = julian_sec - 86400.0 + 2.0
00116   else if ( julian_sec <= 2.0 ) then
00117      julian_dayb(1) = julian_day - 1.0
00118      julian_dayb(2) = julian_day
00119      julian_secb(1) = julian_sec - 2.0
00120      julian_secb(2) = julian_sec + 2.0
00121   else
00122      julian_dayb(1) = julian_day
00123      julian_dayb(2) = julian_day
00124 !rv,sgi
00125 !     julian_secb(1) = julian_sec + 86400.0 - 2.0
00126      julian_secb(1) = julian_sec - 2.0
00127 !rv,sgi
00128      julian_secb(2) = julian_sec + 2.0
00129   endif
00130 
00131 !-----------------------------------------------------------------------
00132 ! 2nd Check for each field in each component whether there are restarts
00133 !     transients to deal with
00134 !-----------------------------------------------------------------------
00135 
00136   do comp_id = 1, Number_of_Comps_allocated
00137 
00138     if ( Comps(comp_id)%status /= PSMILe_status_defined ) cycle
00139 
00140     sga_smioc_transi => sga_smioc_comp(comp_id)%sga_smioc_transi
00141 
00142     nbr_fields_of_comp = 0
00143 
00144     do field_id = 1, Number_of_Fields_allocated
00145        if ( Fields(field_id)%comp_id == comp_id .and. &
00146             Fields(field_id)%status  == PSMILe_status_defined ) &
00147             nbr_fields_of_comp = nbr_fields_of_comp + 1
00148     enddo
00149 
00150     do field_id = 1, nbr_fields_of_comp
00151 
00152        if ( Fields(field_id)%status /= PSMILe_status_defined .or. &
00153             field_id > Number_of_Fields_allocated ) cycle
00154  
00155        fp => Fields(field_id)
00156        smioc_loc     = fp%smioc_loc
00157 
00158        dataType      = fp%dataType
00159        len           = fp%size
00160 
00161        if (smioc_loc == PRISM_UNDEFINED) cycle
00162        if ( .not. associated(sga_smioc_transi(smioc_loc)%sga_transi_out)) cycle
00163 
00164        nb_transiouts = size (sga_smioc_transi(smioc_loc)%sga_transi_out)
00165 
00166        if ( fp%transi_type == PSMILe_bundle ) then
00167           nbr_fields = fp%var_shape(2,Grids(Methods(fp%method_id)%grid_id)%n_dim+1)
00168        else
00169           nbr_fields = 1
00170        endif
00171 
00172 !-----------------------------------------------------------------------
00173 ! 3rd Check for each task within a transient
00174 !-----------------------------------------------------------------------
00175 
00176        do task_id = 1, nb_transiouts
00177 
00178 !-----------------------------------------------------------------------
00179 ! 4th Control lag and existance of restart
00180 !-----------------------------------------------------------------------
00181 !
00182 !      Shift task id into the range of task ids for reading the restart!
00183 !      (Shifted task id is needed by  psmile_read_[int,real,dble])
00184 !      We only enter the section of sending restart data when a lag was set.
00185 !      In this case psmile_check_restart returns 0 or 1 depending on whether
00186 !      a restart was found or not. In case that a lag was defined but no
00187 !      restart is available we send zeros. Perhaps this should only be
00188 !      allowed at the beginning of an experiment for the first job.
00189 !
00190           io_taskid = 2*nb_transiouts+1+task_id
00191 
00192           lag = sga_smioc_transi(smioc_loc)%sga_transi_out(task_id)%ig_lag
00193           if ( lag /= PSMILe_undef .and. lag > 0 ) then
00194              call psmile_check_restart(field_id, task_id, info, ierror)
00195           else
00196              info = -1
00197           endif
00198 
00199           if ( info > 0 ) then
00200 #ifdef VERBOSE
00201              print *, trim(ch_id), ': Reading restarts for field ', trim(fp%local_name)
00202 #endif
00203 
00204 !-----------------------------------------------------------------------
00205 ! 5th Handling of integer fields
00206 !-----------------------------------------------------------------------
00207 
00208              select case (datatype)
00209 
00210              case ( PRISM_Integer )
00211 !
00212 ! a) Check and allocate memory for integers
00213 !
00214                 if ( .not. allocated(idata_array) ) &
00215                      allocate ( idata_array(len), STAT = ierror)
00216                 if ( ierror > 0 ) then
00217                    ierrp (1) = ierror
00218                    ierrp (2) = len
00219 
00220                    ierror = PRISM_Error_Alloc
00221                    call psmile_error ( ierror, 'idata_array', &
00222                         ierrp, 2, __FILE__, __LINE__ )
00223                    return
00224                 endif
00225 !
00226 ! b) Read data from restart
00227 !
00228                 call psmile_read_byid_int ( field_id, io_taskid, idata_array, &
00229                      julian_day, julian_sec, &
00230                      julian_dayb,julian_secb, local_timeop, ierror )
00231 
00232 !
00233 ! c) Send data away
00234 !
00235                 call psmile_put_field_int (field_id, task_id, idata_array, &
00236                      len, nbr_fields, ierror)
00237 
00238 !-----------------------------------------------------------------------
00239 ! 6th Handling of real fields
00240 !-----------------------------------------------------------------------
00241 
00242              case ( PRISM_Real )
00243 !
00244 ! a) Check and allocate memory for reals
00245 !
00246                 if ( .not. allocated(rdata_array) ) &
00247                      allocate ( rdata_array(len), STAT = ierror)
00248                 if ( ierror > 0 ) then
00249                    ierrp (1) = ierror
00250                    ierrp (2) = len
00251 
00252                    ierror = PRISM_Error_Alloc
00253                    call psmile_error ( ierror, 'rdata_array', &
00254                         ierrp, 2, __FILE__, __LINE__ )
00255                    return
00256                 endif
00257 !
00258 ! b) Read data from restart
00259 !
00260                 call psmile_read_byid_real ( field_id, io_taskid, rdata_array, &
00261                      julian_day, julian_sec, &
00262                      julian_dayb,julian_secb, local_timeop, ierror )
00263 !
00264 ! c) Send data away
00265 !
00266                 call psmile_put_field_real (field_id, task_id, rdata_array, &
00267                      len, nbr_fields, ierror)
00268 
00269 !-----------------------------------------------------------------------
00270 ! 7th Handling of double precision fields
00271 !-----------------------------------------------------------------------
00272              
00273              case ( PRISM_Double_Precision )
00274 !
00275 ! a) Check and allocate memory for doubles
00276 !
00277                 if ( .not. allocated(ddata_array) ) &
00278                      allocate ( ddata_array(len), STAT = ierror)
00279                 if ( ierror > 0 ) then
00280                    ierrp (1) = ierror
00281                    ierrp (2) = len
00282 
00283                    ierror = PRISM_Error_Alloc
00284                    call psmile_error ( ierror, 'qdata_array', &
00285                         ierrp, 2, __FILE__, __LINE__ )
00286                    return
00287                 endif
00288 !
00289 ! b) Read data from restart
00290 !
00291                 call psmile_read_byid_dble ( field_id, io_taskid, ddata_array, &
00292                      julian_day, julian_sec, &
00293                      julian_dayb,julian_secb, local_timeop, ierror )
00294 !
00295 ! c) Send data away
00296 !
00297                 call psmile_put_field_dble (field_id, task_id, ddata_array, &
00298                      len, nbr_fields, ierror)
00299 
00300              case default
00301 
00302                 ierror = PRISM_Error_Internal
00303                 call psmile_error ( ierror, 'datatype is currently not supported', &
00304                  ierrp, 0, __FILE__, __LINE__ )
00305 
00306              end select
00307 
00308           else if ( info == 0 ) then
00309 
00310              print *, trim(ch_id), ': WARNING: No restarts for field ', trim(fp%local_name)
00311              call psmile_flushstd()
00312              call psmile_abort
00313 
00314           endif ! info
00315 
00316        enddo ! task_id
00317 
00318     enddo ! field_id
00319 
00320   enddo ! comp_id
00321 
00322 !-----------------------------------------------------------------------
00323 ! 8th Epilogue
00324 !-----------------------------------------------------------------------
00325 
00326   if ( allocated(ddata_array) ) deallocate(ddata_array) 
00327   if ( allocated(rdata_array) ) deallocate(rdata_array) 
00328   if ( allocated(idata_array) ) deallocate(idata_array) 
00329 
00330 #ifdef VERBOSE
00331   print 9990, trim(ch_id), ierror
00332   call psmile_flushstd
00333 
00334 9980 format (1x, a, ': psmile_get_restart: start' )
00335 9990 format (1x, a, ': psmile_get_restart: eof ierror =', i3 )
00336 
00337 #endif /* VERBOSE */
00338 
00339 end subroutine psmile_get_restart

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1