prism_get.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, CERFACS, Toulouse, France.
00003 ! Copyright 2006-2010, SGI Germany, Munich, Germany.
00004 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00005 ! All rights reserved. Use is subject to OASIS4 license terms.
00006 !-----------------------------------------------------------------------
00007 !BOP
00008 !
00009 ! !ROUTINE: PRISM_Get
00010 !
00011 ! !INTERFACE:
00012 
00013    subroutine prism_get ( field_id, date, date_bounds, data_array, info, ierror )
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018       use PRISM_calendar
00019 
00020       use PSMILe
00021       use PSMILe_SMIOC, only : sga_smioc_comp, transient
00022 
00023       Implicit none
00024 !
00025 ! !INPUT PARAMETERS:
00026 !
00027       Integer, Intent (InOut)              :: field_id
00028 
00029 !     Handle to the variable information
00030 
00031       Type(PRISM_Time_Struct), Intent (In) :: date
00032 
00033 !     Date on which the information is located in time
00034 
00035       Type(PRISM_Time_Struct), Intent (In) :: date_bounds(2)
00036 
00037 !     Time interval for which the data is representative
00038 !     lower bound: date_bounds(1), upper bound: date_bounds(2)
00039 !
00040 ! !INPUT/OUTPUT PARAMETERS:
00041 !
00042       Double Precision, Intent (InOut)     :: data_array(*)
00043 
00044 !     the data itself
00045 !     May also contain data which is not changed;
00046 !     therefore Intent (InOut)
00047  
00048 ! !OUTPUT PARAMETERS:
00049 !
00050       Integer, Intent (Out)                :: info
00051 
00052 !     returned info about action performed
00053 
00054       Integer, Intent (Out)                :: ierror
00055 
00056 !     Returns the error code of prism_get;
00057 !             ierror = 0 : No error
00058 !             ierror > 0 : Severe error
00059 !
00060 ! !LOCAL VARIABLES
00061 !
00062       DOUBLE precision            :: my_dble
00063       Double Precision, Parameter :: dble_huge = huge(my_dble)
00064 
00065       Real (PSMILe_float_kind) :: julian_day, julian_dayb(2)
00066       Real (PSMILe_float_kind) :: julian_sec, julian_secb(2)
00067 
00068       Logical            :: action(3)
00069 ! WARNING: WORKAROUND
00070       Logical            :: precise = .false.
00071 
00072       Integer, Parameter :: nerrp = 2
00073       Integer            :: ierrp (nerrp)
00074 
00075       Type (transient), pointer :: sga_smioc_transi(:)
00076 !
00077 ! User-defined interpolation variables
00078 !
00079       Logical            :: ll_userdef
00080 
00081       Type (Userdef), pointer :: ug
00082       Type (GridFunction), pointer :: fp
00083 
00084       Integer            :: il_udef_id
00085       Integer            :: il_dim1, il_dim2, il_nbfld
00086       Integer            :: field_id_1, field_id_2, il_fsize, il_gsppp
00087       Integer            :: il_side
00088       Integer            :: nbr_in
00089       Integer            :: il_size1
00090 !
00091 ! !DESCRIPTION:
00092 !
00093 ! Subroutine "prism_get" get the data from a remote application
00094 ! or the io library. This subroutine performs a blocking get of a
00095 ! variable the PSMILe.
00096 !
00097 !
00098 ! !REVISION HISTORY:
00099 !   Date      Programmer   Description
00100 ! ----------  ----------   -----------
00101 ! 03.07.02    R. Redler    created
00102 !
00103 !EOP
00104 !----------------------------------------------------------------------
00105 !
00106 ! $Id: prism_get.F90 2805 2010-12-07 10:16:28Z coquart $
00107 ! $Author: coquart $
00108 !
00109   Character(len=len_cvs_string), save :: mycvs = 
00110       '$Id: prism_get.F90 2805 2010-12-07 10:16:28Z coquart $'
00111 !
00112 !----------------------------------------------------------------------
00113 
00114 
00115 #ifdef VERBOSE
00116       print *, trim(ch_id), ': prism_get: field_id', field_id
00117 
00118       call psmile_flushstd
00119 #endif /* VERBOSE */
00120 
00121 !-----------------------------------------------------------------------
00122 !  0th Initialization
00123 !-----------------------------------------------------------------------
00124 
00125       ierror = 0
00126       info   = PRISM_NOACTION
00127       ll_userdef = .false.
00128 
00129 !-----------------------------------------------------------------------
00130 ! 1st Preliminary check of field_id and setting the pointer
00131 !-----------------------------------------------------------------------
00132 
00133       if ( field_id == PRISM_UNDEFINED ) then
00134 #ifdef VERBOSE
00135          PRINT*, TRIM(ch_id), ': prism_put we return because field_id undefined:', field_id
00136          print *, trim(ch_id), ': prism_put: eof ierror', ierror
00137          call psmile_flushstd
00138 #endif /* VERBOSE */
00139          return
00140       endif
00141 
00142       if ( field_id > size (Fields) .or. field_id < 1) then
00143 
00144             ierrp (1) = field_id
00145             ierror = PRISM_Error_Arg
00146             call psmile_error ( ierror, 'field_id', &
00147                                 ierrp, 1, __FILE__, __LINE__ )
00148             return
00149       endif
00150 !-----------------------------------------------------------------------
00151 !  2nd Check field_id
00152 !-----------------------------------------------------------------------
00153 
00154       if ( Fields(field_id)%status /= PSMILe_status_defined ) then
00155          ierrp (1) = field_id
00156       
00157          ierror = PRISM_Error_Arg
00158       
00159          call psmile_error ( ierror, 'field_id', &
00160                              ierrp(1), 1, __FILE__, __LINE__ )
00161          return
00162       endif
00163 
00164 !-----------------------------------------------------------------------
00165 ! 3rd Return in case there is nothing to do
00166 !-----------------------------------------------------------------------
00167 
00168       if ( Fields(field_id)%smioc_loc == PRISM_UNDEFINED .or. &
00169          (.not. Fields(field_id)%used_for_coupling      .and. &
00170           .not. Fields(field_id)%used_for_io )) then
00171 #ifdef VERBOSE
00172          print *, trim(ch_id), ': prism_get: nothing to do'
00173          print *, trim(ch_id), ': prism_get: eof ierror ', ierror
00174          call psmile_flushstd
00175 #endif /* VERBOSE */
00176          return
00177       endif
00178 
00179 !-----------------------------------------------------------------------
00180 !  4th check date information
00181 !-----------------------------------------------------------------------
00182 !
00183 !  ... convert date and date bounds into julian days and seconds
00184 !
00185       call psmile_date2ju ( date,           julian_day,     julian_sec    )
00186       call psmile_date2ju ( date_bounds(1), julian_dayb(1), julian_secb(1))
00187       call psmile_date2ju ( date_bounds(2), julian_dayb(2), julian_secb(2))
00188 
00189 !   ... check whether bounds are consistent
00190 
00191       if ( julian_dayb(2) <  julian_dayb(1) .or. &
00192          ( julian_dayb(1) == julian_dayb(2) .and. &
00193            julian_secb(2) <  julian_secb(1) ) ) then
00194 
00195            ierrp (1) = field_id
00196 
00197            ierror = PRISM_Error_Date
00198 
00199            call psmile_error ( ierror, 'upper bound < lower bound', &
00200                              ierrp(1), 1, __FILE__, __LINE__ )
00201            return
00202 
00203       endif
00204 
00205 !  ... check whether date is within bounds
00206 
00207       if ( ( julian_dayb(1) >  julian_day       .or.  &
00208              julian_day     >  julian_dayb(2) ) .or.  &
00209            ( julian_dayb(1) == julian_day       .and. &
00210              julian_sec     <  julian_secb(1) ) .or.  &
00211            ( julian_dayb(2) == julian_day       .and. &
00212              julian_sec     >  julian_secb(2) ) ) then
00213 
00214          ierrp (1) = field_id
00215 
00216          ierror = PRISM_Error_Date
00217 
00218          call psmile_error ( ierror, 'date out of bounds', &
00219                              ierrp(1), 1, __FILE__, __LINE__ )
00220          return
00221 
00222       endif
00223 
00224 !-----------------------------------------------------------------------
00225 !  Preliminary check : user-defined interpolation ?
00226 !-----------------------------------------------------------------------
00227 !
00228       field_id_1 = field_id
00229       fp => Fields(field_id)
00230 
00231 !  Future loop on In_channels
00232       nbr_in = fp%Taskin%nbr_inchannels
00233 #ifdef DEBUG
00234       print *, trim(ch_id),': prism_get: Field_id = ',field_id, nbr_in
00235 #endif /* DEBUG */
00236 #ifdef PRISM_ASSERTION
00237       if (nbr_in == 0) then
00238          print *, trim(ch_id), "prism_get: nbr_in", nbr_in
00239          call PSMILE_Assert (__FILE__, __LINE__, &
00240               "nbr_in == 0")
00241       endif
00242 #endif
00243 
00244 !     do il_i = 1, nbr_in
00245 !        il_udef_id = fp%Taskin%In_channel(il_i)%userdef_id
00246 
00247 !  Current state of the code ... only one channel here
00248       il_udef_id = fp%Taskin%In_channel(1)%userdef_id
00249 #ifdef DEBUG
00250       print *, trim(ch_id),': prism_get: il_udef_id = ',il_udef_id
00251 #endif /* DEBUG */
00252 
00253       if ( il_udef_id /= PSMILe_undef ) then
00254 !
00255 !        data_array received will contain the gridless function
00256          ug => Userdefs(il_udef_id)
00257          ll_userdef = .true.
00258          il_side = ug%ig_transi_side
00259          il_dim1 = size ( fp%var_shape(:,:), dim=1 )
00260          il_dim2 = size ( fp%var_shape(:,:), dim=2 )
00261          field_id_2 = fp%Taskin%In_channel(1)%assoc_var_id
00262 #ifdef DEBUG
00263  PRINT *, ' Field_id, Userdef_id, field_id_2 = ', field_id_1, il_udef_id, field_id_2
00264  PRINT *, ' il_side, number of dims of geo function = ',il_side, il_dim1, il_dim2
00265  PRINT *, ' content of var_shape array = ', fp%var_shape(:,:)
00266 #endif
00267 !        Check length of data
00268 !
00269          il_fsize = Fields(field_id_2)%size       ! computed from actual_shape_pr
00270          il_gsppp = ug%ig_nb_ppp
00271          il_nbfld = ug%ig_nbr_fields              ! defined for prism_def_var
00272 ! Size of a single field (ig_nbr_fields is 1 or nb_bundles)
00273          il_size1 = il_fsize / il_nbfld     ! In case of bundle : size of 1 field
00274 #ifdef DEBUG
00275          PRINT *, ' prism_get : il_nbfld il_fsize = ', il_nbfld, il_fsize
00276          PRINT *, ' prism_get : il_size1 il_gsppp = ', il_size1, il_gsppp
00277 #endif
00278 #ifdef PRISM_ASSERTION
00279       if ( il_size1 /= il_gsppp ) then
00280          call PSMILe_Assert (__FILE__, __LINE__, "Incorrect data size")
00281       endif
00282 #endif
00283 
00284 !        Allocate space for the gridless function
00285 !
00286          if ( Fields(field_id)%dataType == PRISM_Real ) then
00287             Allocate ( ug%real_gridless(1:il_gsppp,1,1,il_nbfld), STAT=ierror )
00288          elseif ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00289             Allocate ( ug%dble_gridless(1:il_gsppp,1,1,il_nbfld), STAT=ierror )
00290          endif
00291 !
00292 !    user-defined interpolation : substitute gridless function and gridless grid
00293 !    From now on : use the gridless function and the gridless grid are in place of
00294 !                  the geographical function and geographical grid
00295 !
00296          Nullify (fp)
00297          field_id = field_id_2
00298 #ifdef DEBUG
00299          PRINT *, ' New (gridless function) field_id = ',field_id
00300 #endif
00301       endif    !  end of code specific to User-defined interpolation
00302 
00303 !  ... check for time gaps between two calls to prism_get
00304 
00305       if ( Fields(field_id)%Taskin%Judate_Ubnd%days /= dble_huge ) then
00306       if ( julian_dayb(1) /= Fields(field_id)%Taskin%Judate_Ubnd%days .or. &
00307            julian_secb(1) /= Fields(field_id)%Taskin%Judate_Ubnd%secs ) then
00308            ierrp (1) = field_id
00309            ierror = PRISM_Error_Date
00310 
00311            print *, trim(ch_id), ': prism_get upper bound from previous get:'
00312            print *, trim(ch_id), ': - days: ', Fields(field_id)%Taskin%Judate_Ubnd%days
00313            print *, trim(ch_id), ': - secs: ', Fields(field_id)%Taskin%Judate_Ubnd%secs
00314 
00315            print *, trim(ch_id), ': prism_get lower bound from this call:'
00316            print *, trim(ch_id), ': - days: ', julian_dayb(1)
00317            print *, trim(ch_id), ': - secs: ', julian_secb(1)
00318 
00319            call psmile_error ( ierror, 'inconsistent date bounds for prism_get', &
00320                                ierrp(1), 1, __FILE__, __LINE__ )
00321            return
00322       endif
00323       endif
00324 
00325 !  ... reset upper and lower bounds 
00326 
00327       Fields(field_id)%Taskin%Judate_Lbnd%days=julian_dayb(1)
00328       Fields(field_id)%Taskin%Judate_Lbnd%secs=julian_secb(1)
00329 
00330       Fields(field_id)%Taskin%Judate_Ubnd%days=julian_dayb(2)
00331       Fields(field_id)%Taskin%Judate_Ubnd%secs=julian_secb(2)
00332 
00333 !------------------------------------------------------------------------
00334 !    5th Check required action
00335 !------------------------------------------------------------------------
00336 
00337       call psmile_check_action ( field_id, 0, precise,       &
00338                                  julian_day, julian_dayb(1), &
00339                                  julian_sec, julian_secb(1), &
00340                                  action )
00341 
00342       if ( .not. action(1) .and. .not. action(2) ) then
00343 #ifdef VERBOSE
00344          print *, trim(ch_id), ': prism_get: eof nothing to do! ierror ', ierror
00345 #endif /* VERBOSE */
00346          return
00347       endif
00348 
00349 !------------------------------------------------------------------------
00350 !    6th Now start getting the data
00351 !------------------------------------------------------------------------
00352 
00353 !     Test if we have a "user-defined" interpolation
00354 
00355       if ( ll_userdef ) then
00356 !  1.    Get the gridless function in gridless structures
00357 !  2.    Restore the geographical field in data_array
00358 
00359          if ( Fields(field_id)%dataType == PRISM_Real ) then
00360             call psmile_get_real ( field_id, julian_day, julian_sec,  &
00361                  julian_dayb, julian_secb, ug%real_gridless, action(1), action(2), &
00362                  info, ierror )
00363             call psmile_gridless_func_real ( field_id_1, il_udef_id, il_side, &
00364                                              data_array, ierror )
00365          endif
00366 
00367          if ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00368             call psmile_get_dble ( field_id, julian_day, julian_sec,  &
00369                  julian_dayb, julian_secb, ug%dble_gridless, action(1), action(2), &
00370                  info, ierror )
00371             call psmile_gridless_func_dble ( field_id_1, il_udef_id, il_side, &
00372                                              data_array, ierror)
00373          endif
00374 
00375 !        Deallocate all memory space of Userdefs structure :  il_udef_id
00376 
00377            if ( Fields(field_id)%dataType == PRISM_Real ) then
00378               Deallocate (ug%real_gridless, STAT=ierror)
00379            else if ( Fields(field_id)%dataType == PRISM_Double_Precision ) then
00380               Deallocate (ug%dble_gridless, STAT=ierror)
00381            endif
00382 
00383            if (ierror > 0) then
00384               ierrp (1) = ierror
00385               ierrp (2) = il_udef_id
00386               ierror = PRISM_Error_Dealloc
00387               call psmile_error ( ierror, 'Userdefs', &
00388                    ierrp, 2, __FILE__, __LINE__ )
00389               return
00390            endif
00391 
00392            Nullify (ug%real_gridless)
00393            Nullify (ug%dble_gridless)
00394 
00395       else
00396 
00397 !  "Normal" case : data_array received IS the geographical function
00398 
00399          if ( Fields(field_id)%dataType == PRISM_Integer ) &
00400             call psmile_get_int ( field_id, julian_day, julian_sec,  &
00401                  julian_dayb, julian_secb, data_array, action(1), action(2), &
00402                  info, ierror )
00403 
00404          if ( Fields(field_id)%dataType == PRISM_Real ) &
00405             call psmile_get_real ( field_id, julian_day, julian_sec,  &
00406                  julian_dayb, julian_secb, data_array, action(1), action(2), &
00407                  info, ierror )
00408 
00409          if ( Fields(field_id)%dataType == PRISM_Double_Precision ) &
00410             call psmile_get_dble ( field_id, julian_day, julian_sec,  &
00411                  julian_dayb, julian_secb, data_array, action(1), action(2), &
00412                  info, ierror )
00413 
00414       endif     ! end of "user-defined" interpolation test
00415 !
00416 !   enddo !    on nbr_in  input channels....
00417 !
00418 
00419 
00420      ! We cannot distinguish here whether data come from a restart or coupler.
00421      if ( action(1) ) info = info + 1000 ! Coupling
00422      if ( action(2) ) info = info +  100 ! IO
00423 
00424 !------------------------------------------------------------------------
00425 !    7th Epilogue
00426 !------------------------------------------------------------------------
00427 
00428 !    Restore input value of field_id  as the geographical field Id
00429 
00430      Nullify(fp)
00431      field_id = field_id_1
00432 
00433 #ifdef VERBOSE
00434       print *, trim(ch_id), ': prism_get: eof ierror ', ierror
00435 
00436       call psmile_flushstd
00437 #endif /* VERBOSE */
00438 !
00439    end subroutine prism_get

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1