psmile_read_byid_int.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, SGI Germany, Munich, Germany.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE:  psmile_read_byid_int
00008 !
00009 ! !INTERFACE:
00010       subroutine psmile_read_byid_int(id_varid,id_taskid,dd_a,ju_day,ju_sec &
00011                                       ,ju_dayb,ju_secb,timeop,ierror)
00012 !
00013 ! !USES:
00014 !
00015       use prism_constants
00016       use psmile, dummy_interface=> psmile_read_byid_int
00017       use prism_calendar
00018 !
00019       implicit none
00020 !
00021 ! !INPUT PARAMETERS:
00022 !
00023       integer,intent(in)::id_varid
00024       integer,intent(in)::id_taskid
00025       integer,target, Intent (inout) :: dd_a(*)
00026       double precision, Intent (In) :: ju_day, ju_sec
00027       double precision, Intent (In) :: ju_dayb(2), ju_secb(2)
00028 !
00029 ! !OUTPUT PARAMETERS:
00030 !
00031       logical,intent(out):: timeop
00032       integer,intent(out):: ierror
00033 !
00034 ! !LOCAL VARIABLES:
00035 !
00036       integer                           :: ierrp(2)
00037       integer               :: il_grid_id,il_grid_type
00038       integer               :: il_method_id,il_varid
00039       integer                           :: il_blockid
00040       integer                           :: iicomp_id
00041       logical                           :: is_block
00042       integer               :: nvcomp,il_ndim,offset,vectorfield
00043       integer                           :: fullsize
00044       integer               :: il_i,lenb,len,ii,jj,il_loc 
00045                                           ,il_time,il_times(2),ibl,ibu
00046       integer               :: myvarshape(2,ndim_3d+2)
00047       integer               :: mygrdshape(2,ndim_3d+2)
00048       integer,allocatable           :: afield(:)
00049       integer                           :: fill
00050       double precision          :: dl_time,dl_ju_day,dl_ju_sec
00051       double precision          :: dl_timeb(2),w,days,sec
00052       double precision,pointer          :: ju_sec_last,ju_day_last 
00053                                          , time_last
00054       double precision,pointer      :: pl_times(:)
00055       double precision              :: aone = 1.d0
00056       logical               :: is_bundle,lmatch,use_domain
00057       logical                           :: lower_in,upper_in
00058       logical,pointer                   :: llast
00059 #ifdef __PSMILE_WITH_IO
00060       Type(IO_cache),pointer            :: pl_cache
00061 #endif
00062       Type(PRISM_Time_struct)           :: cur_date
00063 !
00064 ! !DESCRIPTION:
00065 !  Routine reads a int field variable.
00066 !  It can handle scalar fields and vector fields as well as bundles of
00067 !  the two types of fields.
00068 !  Last but not least it reads subblocks of the above in case of a
00069 !  GME grid.
00070 !  If the current time stamp can not be read by matching the time offset
00071 !  with respect to the current job start date a new file will be opened
00072 !  transparently if the user has given a base name for the input file.
00073 !  The search is based on the extension of the file basename by the job 
00074 !  start dates defining the origin to the time axes contained in those
00075 !  files under consideration.  
00076 !  For a completely specified filename no search is performed.
00077 !  An unsuccessful match of a
00078 !  time stamp leads to a termination of the model component code.
00079 !  
00080 !  The routine honors date bounds for matching time stamps (time lags).
00081 !
00082 !  In case of bundles of vector fields transpose routines are called
00083 !  which transpose the indices of vector components with bundles.
00084 !
00085 ! !REVISION HISTORY:
00086 !
00087 !   Date      Programmer   Description
00088 ! ----------  ----------   -----------
00089 !   21.12.03  R. Vogeslang created
00090 !   29.01.04  R. Vogelsang Time interpolation using date bounds.
00091 !   23.03.04  R. Vogelsang Added the taskid to the interface and related code.
00092 !   12.10.04  R. Vogelsang Values = fill value are set to psmile_undef
00093 !   2.10.05   R. Vogelsang Gaussian reduced grids added
00094 !EOP
00095 !----------------------------------------------------------------------
00096 !
00097       character(len=len_cvs_string),save :: mcvs = 
00098 '$Id: psmile_read_byid_int.F90 2920 2011-01-27 10:55:28Z coquart $'
00099 
00100 !
00101       ierror=0
00102       timeop=.false.
00103 #ifdef __PSMILE_WITH_IO
00104 #ifdef VERBOSE
00105       print*,trim(ch_id),' : PSMILe_read_byid_int: start'
00106       call psmile_flushstd
00107 #endif
00108 !
00109 !     Return from this routine if the field is not subject to file I/O
00110 !     or the I/O action is not MPP_RDONLY (read).
00111 !
00112       if (.not.associated(Fields(id_varid)%io_chan_infos)) return
00113 !
00114 !     IO infos are pointing to io_chan_infos(taskid)
00115 !
00116       if(id_taskid .le. size(Fields(id_varid)%io_task_lookup)) then
00117 
00118         il_i=Fields(id_varid)%io_task_lookup(id_taskid)
00119 
00120       else
00121 
00122         ierror = PRISM_Error_IO_Meta
00123         ierrp (1) = id_taskid
00124         call psmile_error ( ierror &
00125                           , 'Task id out of range! ', &
00126                             ierrp, 1, __FILE__, __LINE__ )
00127         return
00128 
00129       endif
00130 
00131       if (il_i.gt.0) then
00132 
00133          Fields(id_varid)%io_infos =>  Fields(id_varid)%io_chan_infos(il_i)
00134 
00135       else
00136         ierror = PRISM_Error_IO_Meta
00137         ierrp (1) = id_taskid
00138         call psmile_error ( ierror &
00139                           , 'Negative task id! ', &
00140                             ierrp, 1, __FILE__, __LINE__ )
00141         return
00142       endif
00143 
00144       if (Fields(id_varid)%io_infos%action .ne. MPP_RDONLY ) return
00145 !
00146 !     Set the pelist for the current active component
00147 !
00148       iicomp_id=Fields(id_varid)%comp_id
00149       call mpp_set_current_pelist(IO_Comps_infos(iicomp_id)%pelist)
00150 
00151 
00152 
00153       il_varid=id_varid
00154       il_method_id=Fields(id_varid)%method_id
00155       il_grid_id=Methods(il_method_id)%grid_id
00156       il_grid_type=Grids(il_grid_id)%grid_type
00157       il_ndim=Grids(il_grid_id)%n_dim
00158       il_blockid=0
00159       is_block=.false.
00160 
00161       use_domain=.true.
00162       if(Fields(id_varid)%io_infos%threading .eq. MPP_MULTI .and. &
00163          Fields(id_varid)%io_infos%fileset .eq. MPP_MULTI )use_domain=.FALSE.
00164 
00165       dl_ju_day=ju_day
00166       dl_ju_sec=ju_sec
00167 !
00168 !     Time offset with respect to the date of the job start contained
00169 !     in the current file.
00170 !
00171 #ifdef VERBOSE
00172       print*,trim(ch_id),' : PSMILe_read_byid_int: called: ju_start_day= ' &
00173             ,Fields(id_varid)%io_infos%ju_start_day &
00174             ,' ju_start_sec= ',Fields(id_varid)%io_infos%ju_start_sec
00175 #endif
00176       dl_time=86400.d0 *(dl_ju_day - Fields(id_varid)%io_infos%ju_start_day) &
00177              +          (dl_ju_sec - Fields(id_varid)%io_infos%ju_start_sec)
00178 #ifdef VERBOSE
00179       print*,trim(ch_id),' : PSMILe_read_byid_int: called: ju_day= ' &
00180             ,ju_day ,' ju_sec= ',ju_sec,' ju_dayb=',ju_dayb,' ju_secb=' &
00181             ,ju_secb,' offset=',dl_time
00182       call psmile_flushstd
00183 #endif
00184 !
00185 !     Time offset of the date bounds with regards to the job start date
00186 !     read from the currentfile
00187 !
00188       dl_timeb(1:2)=86400.d0 *(ju_dayb(1:2) &
00189                    - Fields(id_varid)%io_infos%ju_start_day) &
00190                    +          (ju_secb(1:2)  &
00191                    - Fields(id_varid)%io_infos%ju_start_sec)
00192 !
00193 !     Search the I/O cache of time stamps for the closed intervall containing
00194 !     the current time offset.
00195 !     
00196       if(associated(Fields(id_varid)%io_infos%p_cache%time_stamps)) then
00197         pl_times=>Fields(id_varid)%io_infos%p_cache%time_stamps
00198         pl_cache=>Fields(id_varid)%io_infos%p_cache
00199       else
00200         ierror=PRISM_Error_IO_Read
00201         call psmile_error(ierror, &
00202               'IO cache for time stamps  not allocated!'  &
00203                         ,ierrp,0, __FILE__, __LINE__)
00204       endif
00205 
00206       lmatch=.FALSE.
00207       llast=>pl_cache%llast
00208       ju_sec_last=>pl_cache%ju_sec_last
00209       ju_day_last=>pl_cache%ju_day_last
00210       time_last=>pl_cache%time_last
00211 !
00212       do while(.not.lmatch)
00213         if(size(pl_times)-1.eq.1) then
00214 !
00215 !       Special treatment if the p_cache contains only 1 time stamp
00216 !
00217           il_i=1
00218           if(abs(dl_time-pl_times(il_i)).lt.1.d-8) then
00219             lmatch=.TRUE.
00220             il_times(1)=il_i
00221             ibl=1
00222             ibu=1
00223             pl_cache%ilast=il_i
00224           else
00225             lower_in=.false.
00226             upper_in=.false.
00227             if(pl_times(il_i).le.dl_timeb(1).and. &
00228                pl_times(il_i).ge.dl_timeb(2)) lower_in=.true.
00229             if(lower_in) then
00230               lmatch=.TRUE.
00231               il_times(1)=il_i
00232               ibl=1
00233               ibu=1
00234               pl_cache%ilast=il_i
00235             else
00236               lmatch=.FALSE.
00237             endif
00238           endif
00239 !
00240         else
00241 !
00242 !       More than 1 time stamp within pcache. Find matching closed intervall.
00243 !
00244           if(llast) then
00245 !
00246 !           From a preceeding file we kept the last time stamp
00247 !
00248             pl_times(0)=86400.d0 *(ju_day_last  &
00249                    - Fields(id_varid)%io_infos%ju_start_day) &
00250                    +          (ju_sec_last  &
00251                    - Fields(id_varid)%io_infos%ju_start_sec)
00252 
00253           else
00254 !
00255 !           Just to skip the time stamp 0 the very first time if no last
00256 !           cached time stamp is present.
00257 !
00258             pl_times(0)=huge(aone)
00259 
00260           endif
00261 
00262 !rv          do il_i=pl_cache%ilast,size(pl_times)-1
00263 !
00264 !         Loop over ubound(time stamps)-1 not  ubound(time stamps)
00265 !
00266           do il_i=pl_cache%ilast,size(pl_times)-2
00267             if ((dl_time.ge.pl_times(il_i)) .and.  &
00268                 (dl_time.le.pl_times(il_i+1)) )exit
00269           enddo
00270 !
00271 !         We don't need the stored data of the preceeding file.
00272 !
00273           if(il_i.gt.1) then
00274             llast=.false.
00275             if(associated(pl_cache%data_int)) &
00276               deallocate(pl_cache%data_int,stat=ierror)
00277             if (ierror > 0) then
00278               ierrp (1) = ierror
00279                ierrp (2) = 1
00280                ierror = PRISM_Error_Alloc
00281                call psmile_error ( ierror, 'deallocate pl_cache%data_int', &
00282                             ierrp, 2, __FILE__, __LINE__ )
00283                return
00284             endif
00285           endif
00286 !
00287 !       Now check the found boundaries of the closed intervall 
00288 !       against the date bounds
00289 !
00290           lower_in=.false.
00291           upper_in=.false.
00292           w=-1.
00293           if(pl_times(il_i).ge.dl_timeb(1).and. &
00294              pl_times(il_i).le.dl_timeb(2)) lower_in=.true.
00295           if(il_i.lt.size(pl_times)-1) then
00296             if(pl_times(il_i+1).ge.dl_timeb(1).and. &
00297                pl_times(il_i+1).le.dl_timeb(2)) upper_in=.true.
00298           endif
00299 !
00300           if(il_i.lt.size(pl_times)-1) then
00301             lmatch=.TRUE.
00302 !
00303 !           1st case: None of the intervall boundaries are inside the
00304 !                     date bounds. Return nearest time step.
00305 !
00306             if((.not.lower_in).and.(.not.upper_in)) then
00307 !
00308 !             Which intervall boundary is nearest?
00309 !
00310               if((dl_time-pl_times(il_i)).lt.(pl_times(il_i+1)-dl_time)) then
00311                 il_times(1)=il_i
00312               else
00313                 il_times(1)=il_i+1
00314               endif
00315               ibl=1
00316               ibu=1
00317 !
00318 !           2nd case: Both intervall boundaries are inside the
00319 !                     date bounds. 
00320 !                     Fields(id_varid)%io_infos%ilag_mode = PSMILe_time_nnghbr:
00321 !                     Return nearest neighbour
00322 !                     Fields(id_varid)%io_infos%ilag_mode = PSMILe_time_linear:
00323 !                     Calculate a weight to return
00324 !                     a weighted mean.
00325 !                     DEFAULT:
00326 !                     Calculate a weight to return
00327 !                     a weighted mean.
00328 !
00329             elseif(lower_in.and.upper_in) then
00330 
00331               Select Case(Fields(id_varid)%io_infos%ilag_mode)
00332               Case (PSMILe_time_nnghbr)
00333                 if((dl_time-pl_times(il_i)).lt.(pl_times(il_i+1)-dl_time)) then
00334                   il_times(1)=il_i
00335                 else
00336                   il_times(1)=il_i+1
00337                 endif
00338                 ibl=1
00339                 ibu=1
00340                 timeop=.true.
00341               Case DEFAULT
00342                 w=pl_times(il_i+1)-pl_times(il_i)
00343                 w=(dl_time-pl_times(il_i))/w
00344                 il_times(1)=il_i
00345                 il_times(2)=il_i+1
00346                 ibl=1
00347                 ibu=2
00348                 timeop=.true.
00349               End Select
00350 !
00351 !           3rd case: Only one of the intervall boundaries is inside the
00352 !                     date bounds. Return what was found.
00353 !
00354             elseif(lower_in) then
00355               il_times(1)=il_i
00356               ibl=1
00357               ibu=1
00358             elseif(upper_in) then
00359               il_times(1)=il_i+1
00360               ibl=1
00361               ibu=1
00362             endif
00363             pl_cache%ilast=il_i
00364 #ifdef VERBOSE
00365       print*,trim(ch_id),' : PSMILe_read_byid_int: Match: ju_day= ',dl_ju_day &
00366             ,' ju_sec= ',dl_ju_sec,' offsets of date bounds=',dl_timeb(1:2) &
00367             ,' time_levels= ',il_times(ibl:ibu),pl_times(il_times(ibl)) &
00368                              ,pl_times(il_times(ibu)),' weight(upper)=',w
00369       call psmile_flushstd
00370 #endif
00371           else
00372             
00373 !
00374 !        We did not find a closed intervall in the current file.
00375 !        Cache all time information and trigger the caching of the
00376 !        data as well.
00377 !
00378             lmatch=.FALSE.
00379             llast=.TRUE.
00380             il_times(1)=il_i
00381             ibl=1
00382             ibu=1
00383             time_last=pl_times(il_i)
00384             ju_sec_last=Fields(id_varid)%io_infos%ju_start_sec+time_last
00385             days=dble(int(ju_sec_last/86400.d0))
00386             ju_sec_last=ju_sec_last-days*86400.d0
00387             ju_day_last=Fields(id_varid)%io_infos%ju_start_day + days
00388 #ifdef VERBOSE
00389       print*,trim(ch_id),' : PSMILe_read_byid_int: No match (EOF): ju_day= ' &
00390             ,dl_ju_day &
00391             ,' ju_sec= ',dl_ju_sec,' offsets of date bounds=',dl_timeb(1:2) &
00392             ,' time_levels= ',il_times(ibl:ibu)
00393       call psmile_flushstd
00394 #endif
00395           endif
00396 !
00397         endif
00398 !
00399 !
00400 !     Vector field ?
00401 !
00402       nvcomp=1
00403       vectorfield=0
00404       if(associated(Methods(il_method_id)%vector_pointer)) then
00405         nvcomp=size(Methods(il_method_id)%vector_pointer%array_of_point_ids)
00406         vectorfield=1
00407       endif
00408 
00409 !
00410 !     Set up shapes depending on the grid type.
00411 !
00412       if ( il_grid_type == PRISM_unstructlonlatvrt) then
00413          len = 1
00414          myvarshape(1:2,1)=Fields(id_varid)%var_shape(1:2,1)
00415          myvarshape(1:2,2)=1
00416          myvarshape(1:2,3)=1
00417          mygrdshape(1:2,1)=Grids(il_grid_id)%grid_shape(1:2,1)
00418          mygrdshape(1:2,2)=1
00419          mygrdshape(1:2,3)=1
00420       else if ( il_grid_type == PRISM_unstructlonlat_sigmavrt .or. &
00421             il_grid_type == PRISM_unstructlonlat_regvrt   .or. &
00422             il_grid_type == PRISM_Gaussreduced_regvrt   .or. &
00423             il_grid_type == PRISM_Gaussreduced_sigmavrt   ) then
00424 
00425          len = 2
00426          myvarshape(1:2,1)=Fields(id_varid)%var_shape(1:2,1)
00427          myvarshape(1:2,2)=1
00428          myvarshape(1:2,3)=Fields(id_varid)%var_shape(1:2,3)
00429 
00430          mygrdshape(1:2,1)=Grids(il_grid_id)%grid_shape(1:2,1)
00431          mygrdshape(1:2,2)=1
00432          mygrdshape(1:2,3)=Grids(il_grid_id)%grid_shape(1:2,3)
00433       else
00434          len = il_ndim
00435          myvarshape(1:2,1:len)=Fields(id_varid)%var_shape(1:2,1:len)
00436          mygrdshape(1:2,1:len)=Grids(il_grid_id)%grid_shape(1:2,1:len)
00437       endif
00438 !
00439 !     Is it a bundle ?
00440 !
00441        
00442       lenb=Fields(id_varid)%var_shape(2,len+vectorfield+1) &
00443          -Fields(id_varid)%var_shape(1,len+vectorfield+1)+1
00444 
00445       if(lenb.gt.1) then
00446         is_bundle=.true.
00447         myvarshape(1:2,4)=Fields(id_varid)%var_shape(1:2,len+vectorfield+1)
00448         mygrdshape(1:2,4)=Fields(id_varid)%var_shape(1:2,len+vectorfield+1)
00449       else ! No bundle
00450         is_bundle=.false.
00451       endif
00452 !
00453 !     Calculate  offset to address vector components.
00454 !     Components are written as separate NetCDF variables with 3 axis
00455 !     declared (Sophie Valckes Nov. 2003)
00456 !
00457       offset=1
00458       do ii = 1, len
00459          offset = offset * (Fields(id_varid)%var_shape(2,ii) &
00460                 - Fields(id_varid)%var_shape(1,ii)+1)
00461       enddo
00462 !
00463 !     For GME like grids with a roundrobin distribution of fields over
00464 !     blocks I match the given varid with the set of related ids.
00465 !     The position of a match will be used as a block id.
00466 !     The file and field  information will be take from the last
00467 !     varid of the set of related ids because each block of the GME grid
00468 !     is partitioned in the same manner.
00469 !
00470       if(associated(Fields(id_varid)%io_infos%related_ids)) then
00471         fullsize=size(Fields(id_varid)%io_infos%related_ids)
00472         if(fullsize.gt.1) then
00473           do il_i=1,fullsize
00474             if(id_varid.eq. Fields(id_varid)%io_infos%related_ids(il_i)) exit
00475           enddo
00476 !
00477 !         The MPP_IO file descriptor is always taken from the first entry
00478 !         of the sets of related ids.
00479 !
00480           il_varid=Fields(id_varid)%io_infos%related_ids(1)
00481 !
00482 !         Take the block number from the matching entry within the set
00483 !         of related ids. 
00484 !
00485           il_blockid=Fields(il_i)%io_infos%block_id
00486           is_block=.true.
00487         endif
00488       endif 
00489 
00490       if(is_bundle) then
00491         offset=offset*lenb
00492         if(nvcomp.gt.1) then
00493           fullsize=offset*nvcomp
00494           allocate(afield(1:fullsize),stat=ierror)
00495           if (ierror > 0) then
00496             ierrp (1) = ierror
00497             ierrp (2) = 1
00498             ierror = PRISM_Error_Alloc
00499             call psmile_error ( ierror, 'afield', &
00500                          ierrp, 2, __FILE__, __LINE__ )
00501             return
00502           endif
00503 
00504         endif
00505       endif
00506 
00507 !
00508 !     Read it
00509 !
00510       
00511       do jj=ibl,ibu
00512       il_time=il_times(jj)
00513 #ifdef VERBOSE
00514       print*,trim(ch_id),' : PSMILe_read_byid_int: ju_day= ',dl_ju_day &
00515             ,' ju_sec= ',dl_ju_sec,' time_level= ',il_time
00516       call psmile_flushstd
00517 #endif
00518       
00519 !
00520 !    Set the first location in the user provided field dd_a.
00521 !
00522 
00523       il_loc=1
00524       
00525       do il_i=1,nvcomp
00526         if(il_time.eq.0.and.llast) then
00527 !
00528 !        We had a match with the last time_stamp of the preceeding file.
00529 !        Return what was cached.
00530 !
00531          dd_a(il_loc:il_loc+offset-1)=pl_cache &
00532                                      %data_int(il_loc:il_loc+offset-1)
00533         else
00534         
00535         ii=Fields(il_varid)%io_infos%p_mpp_io%findex(il_i)
00536 
00537         if(.not.is_bundle) then
00538         call psmile_read_3d_int(Fields(il_varid)%io_infos%file_unit &
00539                                  ,Fields(il_varid)%io_infos%p_mpp_io%field(ii) &
00540                                  ,Fields(il_varid)%io_infos%p_mpp_io%domain(1) &
00541                                  ,dd_a(il_loc) &
00542                                  ,myvarshape &
00543                                  ,mygrdshape &
00544                                  ,il_time,.TRUE.,il_blockid,is_block &
00545                                  ,use_domain,ierror)
00546         else if (is_bundle) then
00547           if(vectorfield.eq.0) then
00548             call psmile_read_4d_int(Fields(il_varid)%io_infos%file_unit &
00549                                  ,Fields(il_varid)%io_infos%p_mpp_io%field(ii) &
00550                                  ,Fields(il_varid)%io_infos%p_mpp_io%domain(1) &
00551                                      ,dd_a(il_loc) &
00552                                      ,myvarshape &
00553                                      ,mygrdshape &
00554                                      ,il_time,.TRUE.,il_blockid,is_block &
00555                                      ,use_domain,ierror)
00556           else
00557             call psmile_read_4d_int(Fields(il_varid)%io_infos%file_unit &
00558                                  ,Fields(il_varid)%io_infos%p_mpp_io%field(ii) &
00559                                  ,Fields(il_varid)%io_infos%p_mpp_io%domain(1) &
00560                                      ,afield(il_loc) &
00561                                      ,myvarshape &
00562                                      ,mygrdshape &
00563                                      ,il_time,.TRUE.,il_blockid,is_block &
00564                                      ,use_domain,ierror)
00565           endif
00566         endif
00567         endif
00568         il_loc=il_loc+offset
00569       enddo
00570 
00571       if(is_bundle) then
00572         if(nvcomp.gt.1) then
00573 !
00574 !         Individual vector components as bundles were read by looping over
00575 !         vector indices as the last dimension.
00576 !         However, in PSMILE the bundles have to appear as the last 
00577 !         set of indices.
00578 !         Therefore, a transpose of the last two indices has to take place.
00579 !
00580 
00581           if(len.eq.1) call trs_vec_bundle_1d_int(afield,dd_a &
00582                                                   ,Fields(id_varid)%var_shape &
00583                                                   ,mygrdshape &
00584                                                   ,ierror)
00585           if(len.eq.2) call trs_vec_bundle_2d_int(afield,dd_a &
00586                                                   ,Fields(id_varid)%var_shape &
00587                                                   ,mygrdshape &
00588                                                   ,ierror)
00589           if(len.eq.3) call trs_vec_bundle_3d_int(afield,dd_a &
00590                                                   ,Fields(id_varid)%var_shape &
00591                                                   ,mygrdshape &
00592                                                   ,ierror)
00593         endif
00594       endif
00595 !
00596 !
00597 !     In case that the end of the time stamps in one file was hit or
00598 !     or two time stamps fall into the date bounds we have to cache one
00599 !     set of data.
00600 !
00601         if((jj.eq.1.and.ibu.eq.2).or.llast) then
00602 #ifdef DEBUG
00603       print*,trim(ch_id),' : PSMILe_read_byid_int:',nvcomp*offset &
00604             ,' integers must be cached . Allocating pcache!'
00605       call psmile_flushstd
00606 #endif
00607 !
00608           if(associated(pl_cache%data_int)) then
00609              if(size(pl_cache%data_int).lt.nvcomp*offset) &
00610                 deallocate(pl_cache%data_int,stat=ierror)
00611              if (ierror > 0) then
00612                ierrp (1) = ierror
00613                ierrp (2) = 1
00614                ierror = PRISM_Error_Alloc
00615                call psmile_error ( ierror, 'deallocate pl_cache%data_int', &
00616                              ierrp, 2, __FILE__, __LINE__ )
00617                return
00618              endif
00619           endif
00620 
00621 !
00622           if(.not.associated(pl_cache%data_int)) &
00623              allocate(pl_cache%data_int(1:nvcomp*offset),stat=ierror)
00624           if (ierror > 0) then
00625               ierrp (1) = ierror
00626               ierrp (2) = 1
00627               ierror = PRISM_Error_Alloc
00628               call psmile_error ( ierror, 'allocate pl_cache%data_int', &
00629                             ierrp, 2, __FILE__, __LINE__ )
00630               return
00631           endif
00632 !
00633           pl_cache%data_int(1:nvcomp*offset)=dd_a(1:nvcomp*offset)
00634 #ifdef DEBUG
00635       print*,trim(ch_id),' : PSMILe_read_byid_int:',nvcomp*offset &
00636             ,' ints cached'
00637       call psmile_flushstd
00638 #endif
00639         endif
00640 !
00641       enddo  ! End of loop over matching time stamps
00642 !
00643 !     Compute the waited mean of the time stamps i and i+1
00644 !
00645       if(ibu.eq.2.and.lmatch) then
00646         dd_a(1:nvcomp*offset)=nint((aone-w)*pl_cache%data_int(1:nvcomp*offset) &
00647                              +w*dd_a(1:nvcomp*offset))
00648       endif
00649 
00650       if (.not.lmatch) then
00651 !
00652 !         Convert Julian date to Gregorian date in order to construct 
00653 !         a file name if a search for a file containing a new, 
00654 !         advanced job start date has
00655 !         to be performed in order to match the  current date.
00656 !
00657           call psmile_ju2date ( cur_date, dl_ju_day, dl_ju_sec )
00658 !   
00659 !         Check file size and open transparently a new file if necessary.
00660 !
00661           call psmile_open_file_byid(id_varid,id_taskid,cur_date,ierror)
00662       endif
00663       enddo   !End of do while
00664 !
00665 !     Save some memory!
00666 !
00667       if(.not.llast) then
00668         if(associated(pl_cache%data_int)) &
00669            deallocate(pl_cache%data_int,stat=ierror)
00670         if (ierror > 0) then
00671           ierrp (1) = ierror
00672           ierrp (2) = 1
00673           ierror = PRISM_Error_Alloc
00674           call psmile_error ( ierror, 'deallocate pl_cache%data_int', &
00675                               ierrp, 2, __FILE__, __LINE__ )
00676           return
00677         endif
00678       endif
00679       
00680 
00681 #ifdef VERBOSE
00682       print*,trim(ch_id),' : PSMILe_read_byid_int: end'
00683       call psmile_flushstd
00684 #endif
00685 #endif
00686 
00687       end subroutine psmile_read_byid_int
00688 
00689       subroutine trs_vec_bundle_1d_int(user,xio,shape,vshape,ierror)
00690       integer,intent(out)         :: ierror
00691       integer,intent(in)      :: shape(2,3)
00692       integer,intent(in)      :: vshape(2,3)
00693       integer,intent(in)          :: user(shape(1,1):shape(2,1) 
00694                                          ,shape(1,3):shape(2,3) 
00695                                          ,shape(1,2):shape(2,2))
00696       integer,intent(out)         ::  xio(shape(1,1):shape(2,1) 
00697                                          ,shape(1,2):shape(2,2) 
00698                                          ,shape(1,3):shape(2,3))
00699 !
00700 !     Local
00701 !
00702       integer                     :: i,j
00703       
00704       ierror=0
00705 
00706       do j=vshape(1,3),vshape(2,3)
00707         do i=vshape(1,2),vshape(2,2)
00708           xio(vshape(1,1):vshape(2,1),i,j)=user(vshape(1,1):vshape(2,1),j,i)
00709         enddo
00710       enddo
00711       
00712       end subroutine trs_vec_bundle_1d_int
00713 
00714       subroutine trs_vec_bundle_2d_int(user,xio,shape,vshape,ierror)
00715       integer,intent(out)         :: ierror
00716       integer,intent(in)      :: shape(2,4)
00717       integer,intent(in)      :: vshape(2,4)
00718       integer,intent(in)          :: user(shape(1,1):shape(2,1) 
00719                                          ,shape(1,2):shape(2,2) 
00720                                          ,shape(1,4):shape(2,4) 
00721                                          ,shape(1,3):shape(2,3))
00722       integer,intent(out)         ::  xio(shape(1,1):shape(2,1) 
00723                                          ,shape(1,2):shape(2,2) 
00724                                          ,shape(1,3):shape(2,3) 
00725                                          ,shape(1,4):shape(2,4))
00726 !
00727 !     Local
00728 !
00729       integer                     :: i,j
00730       
00731       ierror=0
00732 
00733       do j=vshape(1,4),vshape(2,4)
00734         do i=vshape(1,3),vshape(2,3)
00735 
00736           xio(vshape(1,1):vshape(2,1) &
00737              ,vshape(1,2):vshape(2,2),i,j)=user(vshape(1,1):vshape(2,1) &
00738                                                ,vshape(1,2):vshape(2,2), j,i)
00739         enddo
00740       enddo
00741       
00742       end subroutine trs_vec_bundle_2d_int
00743 
00744       subroutine trs_vec_bundle_3d_int(user,xio,shape,vshape,ierror)
00745       integer,intent(out)         :: ierror
00746       integer,intent(in)      :: shape(2,5)
00747       integer,intent(in)      :: vshape(2,5)
00748       integer,intent(in)          :: user(shape(1,1):shape(2,1) 
00749                                          ,shape(1,2):shape(2,2) 
00750                                          ,shape(1,3):shape(2,3) 
00751                                          ,shape(1,5):shape(2,5) 
00752                                          ,shape(1,4):shape(2,4))
00753       integer,intent(out)         ::  xio(shape(1,1):shape(2,1) 
00754                                          ,shape(1,2):shape(2,2) 
00755                                          ,shape(1,3):shape(2,2) 
00756                                          ,shape(1,4):shape(2,4) 
00757                                          ,shape(1,5):shape(2,5))
00758 !
00759 !     Local
00760 !
00761       integer                     :: i,j
00762       
00763       ierror=0
00764 
00765       do j=vshape(1,5),vshape(2,5)
00766         do i=vshape(1,4),vshape(2,4)
00767 
00768           xio(vshape(1,1):vshape(2,1) &
00769              ,vshape(1,2):vshape(2,2) &
00770              ,vshape(1,3):vshape(2,3),i,j)=user(vshape(1,1):vshape(2,1) &
00771                                                ,vshape(1,2):vshape(2,2) &
00772                                                ,vshape(1,3):vshape(2,3) , j,i)
00773         enddo
00774       enddo
00775       
00776       end subroutine trs_vec_bundle_3d_int

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1