psmile_get_int.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_Get_int 
00008 !
00009 ! !INTERFACE:
00010 
00011 subroutine psmile_get_int  ( field_id, julian_day, julian_sec, &
00012               julian_dayb, julian_secb, data_array, cpl_action, io_action, &
00013               info, ierror )
00014 !
00015 ! !USES:
00016 !
00017   use PRISM_constants
00018 
00019   use PSMILe
00020   use PSMILe_SMIOC, only : sga_smioc_comp, transient, transient_in
00021 
00022   implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026   Integer, Intent (In)                 :: field_id
00027 
00028 !     Handle to the variable information
00029 
00030   Real (PSMILe_float_kind), Intent (In)   :: julian_day, julian_dayb(2)
00031   Real (PSMILe_float_kind), Intent (In)   :: julian_sec, julian_secb(2)
00032 
00033 !     Date with  bounds on which the information is located in time
00034 !     Time interval for which the data is representative
00035 !     lower bound: julian_???(1), upper bound: julian_???b(2)
00036 !
00037   Logical, Intent (In)            :: cpl_action
00038 
00039   Logical, Intent (In)            :: io_action
00040 !
00041 !  Switches to activate calls for coupling and/or I/O
00042 !
00043 ! !OUTPUT PARAMETERS:
00044 !
00045   Integer, Intent (InOut)          :: data_array(*)
00046 
00047 !     The data itself
00048 !
00049   Integer, Intent (InOut)          :: info
00050 
00051 !     Returned info about action performed
00052 
00053   Integer, Intent (Out)            :: ierror
00054 
00055 !     Returns the error code of psmile_get_int ;
00056 !             ierror = 0 : No error
00057 !             ierror > 0 : Severe error
00058 !
00059 ! !LOCAL VARIABLES
00060 !
00061   Integer, Allocatable           :: data_scattered(:)
00062 
00063   Integer            :: len_gathered, len_scattered, i, j
00064 
00065   Integer            :: len, len_3d
00066 
00067   Integer            :: nbr_fields ! number of fields in a bundle
00068 
00069   Logical            :: local_timeop
00070   Logical            :: local_add
00071   Logical            :: local_multiply
00072   Logical            :: local_gather
00073 
00074   Integer            :: task_id
00075 
00076   Integer            :: n_grid_dim
00077 
00078   Type (GridFunction), Pointer :: fp
00079 
00080 ! SMIOC data
00081 
00082   Type (transient),    Pointer :: sga_smioc_transi(:)
00083   Type (transient_in), Pointer :: sg_transi_in
00084 
00085 ! For I/O
00086 
00087 #ifdef __PSMILE_WITH_IO
00088   Logical            :: debug_action
00089   Integer            :: il_taskid,il_transiouts,il_smioc_loc
00090 #endif
00091 
00092 ! For Statistics
00093 
00094   Integer            :: stat_nsum
00095   Integer            :: nsum
00096 
00097   Integer            :: shape_in(2,6)
00098   Integer            :: shape_out(2,6)
00099 
00100   Logical            :: mask_set
00101   Logical            :: mask_needed
00102   Logical, Pointer   :: mask_array(:)
00103   Logical, Pointer   :: mask_aux(:)
00104 
00105   Integer, Allocatable          :: recv_buf(:)
00106   Integer, Allocatable          :: stat_max(:)
00107   Integer, Allocatable          :: stat_min(:)
00108   Integer, Allocatable          :: stat_sum(:)
00109   Integer, Allocatable          :: stat_mean(:)
00110 
00111   Character(len=6), Save :: cstats (3)
00112 
00113 ! Error parameters
00114 
00115   Integer, Parameter :: nerrp = 2
00116   Integer            :: ierr, ierrp (nerrp)
00117 !
00118 ! !TO DO:
00119 !
00120 !   Apply error handling for coupling and I/O
00121 !   Apply gathering for bundles of vectors
00122 !
00123 ! !DESCRIPTION:
00124 !
00125 ! Subroutine "psmile_get_int " get the data from a remote application
00126 !                        or the io library
00127 !
00128 !
00129 ! !REVISION HISTORY:
00130 !
00131 !   Date      Programmer   Description
00132 ! ----------  ----------   -----------
00133 ! 03.07.02    R. Redler    created
00134 !
00135 !EOP
00136 !----------------------------------------------------------------------
00137 !
00138 ! $Id: psmile_get_int.F90 2687 2010-10-28 15:15:52Z coquart $
00139 ! $Author: coquart $
00140 !
00141   Character(len=len_cvs_string), save :: mycvs = 
00142        '$Id: psmile_get_int.F90 2687 2010-10-28 15:15:52Z coquart $'
00143 !
00144   Data cstats /'masked', 'valid', 'all'/
00145 !
00146 !----------------------------------------------------------------------
00147 !
00148 #ifdef VERBOSE
00149   print *, trim(ch_id), ': psmile_get_int : field_id', field_id
00150 
00151   call psmile_flushstd
00152 #endif /* VERBOSE */
00153 
00154 !-----------------------------------------------------------------------
00155 !  1st Initialization
00156 !-----------------------------------------------------------------------
00157 
00158    ierror         = 0
00159    task_id        = 1
00160    local_timeop   = .false.
00161 
00162    fp      => Fields(field_id)
00163    sga_smioc_transi => sga_smioc_comp(Fields(field_id)%comp_id)%sga_smioc_transi
00164    sg_transi_in  => sga_smioc_transi(fp%smioc_loc)%sg_transi_in
00165 
00166    n_grid_dim = Grids(Methods(fp%method_id)%grid_id)%n_dim
00167 
00168 !-----------------------------------------------------------------------
00169 ! Determine size of data_array
00170 !-----------------------------------------------------------------------
00171 
00172    len = fp%size
00173 
00174    if ( fp%transi_type == PSMILe_bundle ) then
00175        nbr_fields = fp%var_shape(2,n_grid_dim+1)
00176    else
00177        nbr_fields = 1
00178    endif
00179 
00180    len_3d = len / nbr_fields
00181 
00182 !-----------------------------------------------------------------------
00183 ! Set switches
00184 !-----------------------------------------------------------------------
00185 
00186    local_add      = sg_transi_in%sg_tgt_local_trans%dg_add_scalar  /= PSMILe_dundef
00187    local_multiply = sg_transi_in%sg_tgt_local_trans%dg_mult_scalar /= PSMILe_dundef
00188    local_gather   = sg_transi_in%sg_tgt_local_trans%ig_gather == PSMILe_gath
00189 
00190    mask_set = fp%mask_id /= PRISM_UNDEFINED
00191 
00192    if (mask_set) mask_set = Masks(fp%mask_id)%status /= PSMILe_status_free
00193 
00194 #ifdef __PSMILE_WITH_IO
00195 
00196 !-----------------------------------------------------------------------
00197 !  2nd Perform I/O if required
00198 !-----------------------------------------------------------------------
00199 
00200    if ( io_action ) then
00201 
00202       if ( local_gather ) then
00203 
00204          call psmile_read_byid_int  ( field_id, task_id, data_scattered, &
00205                                       julian_day, julian_sec, &
00206                                       julian_dayb,julian_secb,local_timeop,ierror )
00207       else
00208 
00209          call psmile_read_byid_int  ( field_id, task_id, data_array, &
00210                                       julian_day, julian_sec, &
00211                                       julian_dayb,julian_secb,local_timeop,ierror )
00212       endif
00213 
00214       if ( local_timeop ) info = info + 1
00215 
00216    endif
00217 
00218 #endif
00219 
00220 !-----------------------------------------------------------------------
00221 !  3rd Receive data if required
00222 !-----------------------------------------------------------------------
00223 
00224    if ( cpl_action ) then
00225  
00226       if ( local_gather ) then
00227 
00228          call psmile_get_field_int  (field_id, data_scattered, len_3d, &
00229                                      nbr_fields, ierror)
00230       else
00231 
00232          call psmile_get_field_int  (field_id, data_array, len_3d, &
00233                                      nbr_fields, ierror)
00234       endif
00235 
00236    endif
00237 
00238 !-----------------------------------------------------------------------
00239 !  Perform local transformation after I/O or transfer if required
00240 !-----------------------------------------------------------------------
00241    !
00242    ! ... Multiplication with a constant
00243    !
00244    IF ( local_multiply ) THEN
00245 
00246        IF ( local_gather ) THEN
00247 
00248            data_scattered(1:len) = data_scattered(1:len) * &
00249               sg_transi_in%sg_tgt_local_trans%dg_mult_scalar
00250 
00251        ELSE
00252 
00253           data_array(1:len) = data_array(1:len) * &
00254              sg_transi_in%sg_tgt_local_trans%dg_mult_scalar
00255 
00256       ENDIF
00257 
00258   ENDIF
00259 !  
00260    ! ... Addition of a constant
00261    !
00262    IF ( local_add ) THEN
00263 
00264        IF ( local_gather ) THEN
00265 
00266            data_scattered(1:len) = data_scattered(1:len) + &
00267               sg_transi_in%sg_tgt_local_trans%dg_add_scalar
00268        ELSE
00269 
00270             data_array(1:len) = data_array(1:len) + &
00271                  sg_transi_in%sg_tgt_local_trans%dg_add_scalar
00272        ENDIF
00273        !
00274    ENDIF
00275    !
00276    !
00277 
00278 #ifdef __PSMILE_WITH_IO
00279 !
00280 !  Debug mode
00281 !
00282   il_smioc_loc=fp%smioc_loc
00283   debug_action= sg_transi_in%ig_debugmode .eq. PSMILe_true
00284   if ( debug_action ) then
00285        il_transiouts=0
00286        if(associated(sga_smioc_transi(il_smioc_loc)%sga_transi_out)) &
00287          il_transiouts=size(sga_smioc_transi(il_smioc_loc)%sga_transi_out)
00288 !
00289 !      Shift task id into the range of task ids for debugging!
00290 !
00291        il_taskid=2*il_transiouts+task_id
00292 !
00293 !      Coherence check
00294 !
00295        if(il_taskid.le.size( fp%io_task_lookup)) &
00296          debug_action= fp%io_task_lookup(il_taskid).gt.0
00297   endif
00298 
00299 #ifdef VERBOSE
00300    print *, trim(ch_id), ': psmile_get_int : debug_action', debug_action
00301 
00302    call psmile_flushstd
00303 #endif /* VERBOSE */
00304 
00305 !------------------------------------------------------------------------
00306 ! Writing of debug fields considering local transformations if required
00307 !------------------------------------------------------------------------
00308 !
00309   if ( debug_action ) then
00310 
00311        if ( local_gather ) then
00312 
00313           call psmile_write_byid_int  ( field_id, il_taskid, data_scattered, &
00314                julian_day, julian_sec, ierror )
00315        else
00316 
00317 
00318           call psmile_write_byid_int  ( field_id, il_taskid, data_array, &
00319                 julian_day, julian_sec, ierror )
00320 
00321        endif ! local_gather
00322 
00323    endif ! debug_action
00324 
00325 #endif
00326 
00327 
00328 !-----------------------------------------------------------------------
00329 !  4th Gather data if necessary
00330 !-----------------------------------------------------------------------
00331 
00332    if ( local_gather ) then
00333 
00334 !  ... Get lengths of source and target arrays. The lengths exclude the number of
00335 !      subgrids or physical fields stored in the array. This information is
00336 !      provided with an addition parameter nbr_fields. 
00337 
00338        if ( .not. mask_set ) then
00339           
00340           len_gathered  = len_3d
00341           len_scattered = len_gathered
00342 
00343        else
00344 
00345           len_gathered  = 0
00346           len_scattered = 1
00347 
00348           do i = 1, n_grid_dim
00349 
00350              len_scattered = len_scattered                                &
00351                   * ( Masks(fp%mask_id)%mask_shape(2,i) -   &
00352                   Masks(fp%mask_id)%mask_shape(1,i) + 1 )
00353 
00354              if ( len_scattered /= len_3d ) then
00355 
00356                 ierror = PRISM_Error_Parameter
00357 
00358                 ierrp(1) = len_scattered
00359                 ierrp(2) = len_3d
00360 
00361                 call psmile_error ( PRISM_Error_Parameter, 'Size of Mask and Field', &
00362                      ierrp, nerrp, __FILE__, __LINE__ )
00363 
00364                 return
00365              endif
00366 
00367           enddo
00368 
00369           len_gathered = count ( Masks(fp%mask_id)%mask_array )
00370 
00371        endif
00372 
00373 !  ... now gather scattered data onto the data_array
00374 
00375           call psmile_loc_trans_int  ( PSMILe_gath, nbr_fields, &
00376                                       len_scattered, data_scattered,   &
00377                                       len_gathered, data_array, field_id )
00378    endif
00379 
00380 !-----------------------------------------------------------------------
00381 ! 5th  Calculate statistics
00382 !
00383 !      j = 1: for masked points
00384 !      j = 2: for valid points
00385 !      j = 3: for all points
00386 !
00387 !-----------------------------------------------------------------------
00388 
00389    mask_needed = sg_transi_in%iga_stats(1) == PSMILe_true .or. &
00390                  sg_transi_in%iga_stats(2) == PSMILe_true .or. &
00391                  sg_transi_in%iga_stats(3) == PSMILe_true
00392 
00393    if ( mask_needed ) then
00394 
00395 !     Get dimensions for Statistics vectors
00396       
00397       shape_out      = 1
00398       shape_out(2,6) = nbr_fields
00399 
00400       shape_in = 1
00401 
00402       do i = 1, n_grid_dim
00403          shape_in(1,i) = fp%var_shape(1,i)
00404          shape_in(2,i) = fp%var_shape(2,i)
00405       enddo
00406 
00407       shape_in(2,6)=nbr_fields
00408 
00409 !     Allocate work buffers for statistics
00410 
00411       allocate (mask_aux(len/nbr_fields), STAT=ierr)
00412       if ( ierr /= 0 ) then
00413          ierrp (1) = 1
00414          ierror = PRISM_Error_Alloc
00415          call psmile_error ( ierror, 'mask_aux', ierrp(1), 1, &
00416                              __FILE__, __LINE__ )
00417          return
00418       endif
00419 
00420       allocate (recv_buf(1:shape_out(2,6)), &
00421                 stat_max(1:shape_out(2,6)), &
00422                 stat_min(1:shape_out(2,6)), &
00423                 stat_sum(1:shape_out(2,6)), &
00424                 stat_mean(1:shape_out(2,6)), STAT=ierror)
00425       if ( ierror /= 0 ) then
00426          ierrp (1) = nbr_fields * 5
00427          ierror = PRISM_Error_Alloc
00428          call psmile_error ( ierror, 'recv_buf, stat_{max,min,sum,mean}', &
00429                     ierrp(1), 1, __FILE__, __LINE__ )
00430          return
00431       endif
00432 
00433 !     For all statistics currently supported
00434 
00435       do j = 1, 3
00436 
00437       if ( sg_transi_in%iga_stats(j) == PSMILe_true ) then
00438 
00439          !-----------------------------------------------------------------------
00440          ! In case that there is no mask defined for the field create one.
00441          !-----------------------------------------------------------------------
00442 
00443          if ( mask_set ) then
00444 
00445             select case (j)
00446             case (1)
00447                mask_aux = .not. Masks(fp%mask_id)%mask_array
00448                mask_array => mask_aux            
00449             case (2)
00450                mask_array => Masks(fp%mask_id)%mask_array
00451             case (3)
00452                mask_aux = .true.
00453                mask_array => mask_aux            
00454             end select
00455 
00456          else
00457 
00458             mask_aux = j > 1
00459             mask_array => mask_aux            
00460          endif
00461 
00462          if ( local_gather ) then
00463 
00464             call psmile_multi_reduce_int  ( PSMILe_max, shape_in, data_scattered, &
00465                  shape_out, stat_max, mask_array, ierror )
00466 
00467             call psmile_multi_reduce_int  ( PSMILe_min, shape_in, data_scattered, &
00468                  shape_out, stat_min, mask_array, ierror )
00469 
00470             call psmile_multi_reduce_int  ( PSMILe_integral, shape_in, data_scattered, &
00471                  shape_out, stat_sum, mask_array, ierror )
00472          else
00473 
00474             call psmile_multi_reduce_int  ( PSMILe_max, shape_in, data_array, &
00475                  shape_out, stat_max, mask_array, ierror )
00476 
00477             call psmile_multi_reduce_int  ( PSMILe_min, shape_in, data_array, &
00478                  shape_out, stat_min, mask_array, ierror )
00479 
00480             call psmile_multi_reduce_int  ( PSMILe_integral, shape_in, data_array, &
00481                  shape_out, stat_sum, mask_array, ierror )
00482 
00483          endif
00484 
00485          if ( Comps(fp%comp_id)%act_comm /= MPI_COMM_NULL ) then
00486 
00487             call MPI_Allreduce ( stat_max, recv_buf, shape_out(2,6), MPI_INTEGER, &
00488                  MPI_MAX, Comps(fp%comp_id)%act_comm, ierror )
00489             stat_max (:) = recv_buf (:)
00490 
00491             call MPI_Allreduce ( stat_min, recv_buf, shape_out(2,6), MPI_INTEGER, &
00492                  MPI_MIN, Comps(fp%comp_id)%act_comm, ierror )
00493             stat_min (:) = recv_buf (:)
00494 
00495             call MPI_Allreduce ( stat_sum, recv_buf, shape_out(2,6), MPI_INTEGER, &
00496                  MPI_SUM, Comps(fp%comp_id)%act_comm, ierror )
00497 
00498             stat_sum (:) = recv_buf (:)
00499 
00500             !if the whole field is being used
00501             if (j == 3) then
00502                stat_nsum = len_3d
00503             else
00504                stat_nsum = count(mask_array)
00505             endif
00506 
00507             call MPI_Allreduce ( stat_nsum, nsum, 1, MPI_INTEGER, &
00508                  MPI_SUM, Comps(fp%comp_id)%act_comm, ierror )
00509 
00510             if (nsum > 0) then
00511                stat_mean(:) = stat_sum(:) / nsum
00512             else
00513                stat_mean = 0.0
00514             endif
00515 
00516          endif
00517 
00518 !        Print out statistical data
00519 
00520          write (*, 9990) trim(ch_id)
00521 
00522          write (*, 9980) trim(ch_id), trim(cstats (j))
00523 
00524          write (*,*) trim(ch_id), &
00525               ': ... for field ', trim(fp%local_name)
00526 
00527          write (*,'(1x,a,a)') trim(ch_id), &
00528               ': BundleNr.       Min           Max           Sum           Mean'
00529 
00530          write (*, 9950) trim(ch_id)
00531 
00532          do i = 1, shape_out(2,6)
00533             write(*,'(1x,a,a2,i3,6x,4(1x,e14.6))') trim(ch_id), &
00534                  ': ', i, stat_min(i), stat_max(i), stat_sum(i), stat_mean(i)
00535          enddo
00536 
00537          write (*, 9990) trim(ch_id)
00538 
00539          call psmile_flushstd
00540 
00541       endif ! iga_stats
00542 
00543       enddo ! j-loop
00544 
00545 !  ... Deallocate work array allocated for statistics
00546 
00547       deallocate (recv_buf, stat_min, &
00548                   stat_max, stat_sum, stat_mean, STAT = ierror)
00549 #if defined ( VERBOSE) 
00550       if ( ierror /= 0 ) then
00551          ierrp (1) = nbr_fields
00552          ierror = PRISM_Error_Dealloc
00553          call psmile_error ( ierror, 'recv_buf, stat_{min,max,sum,mean}', &
00554               ierrp(1), 1, __FILE__, __LINE__ )
00555          return
00556       endif
00557 #endif
00558 
00559       deallocate (mask_aux, STAT = ierror)
00560 #if defined ( VERBOSE) 
00561       if ( ierror /= 0 ) then
00562          ierrp (1) = 1
00563          ierror = PRISM_Error_Dealloc
00564          call psmile_error ( ierror, 'mask_aux', &
00565               ierrp(1), 1, __FILE__, __LINE__ )
00566          return
00567       endif
00568 #endif
00569 
00570    endif
00571 
00572 
00573 #ifdef VERBOSE
00574       print *, trim(ch_id), ': psmile_get_int : eof ierror ', ierror
00575 
00576       call psmile_flushstd
00577 #endif /* VERBOSE */
00578 
00579 !  Formats
00580 
00581 9990  format (1x, a, ': ', 65('='))
00582 9980  format (1x, a, ': Statistics over ', a, ' points')
00583 9950  format (1x, a, ': ', 65('-'))
00584 
00585 end subroutine psmile_get_int 

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1