psmile_recv_req_coords_real.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_Recv_req_coords_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_recv_req_coords_real (sender, tag,      &
00012                                               search, recv_req, &
00013                                               receive_mask, new_search, ierror)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use PSMILe, dummy_interface => PSMILe_Recv_req_coords_real
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       Integer, Intent (In)                :: sender
00026 
00027 !     Specifies the sender of the subgrid.
00028 
00029       Integer, Intent (In)                :: tag
00030 
00031 !     Specifies the message tag used
00032 
00033       Logical, Intent (In)                :: receive_mask
00034 
00035 !     Should mask array be received ?
00036 !
00037       Logical, Intent (In)                :: new_search
00038 
00039 !     Option:
00040 !     new_search = .true.  : New search; fill search structure
00041 !     new_search = .fasle. : Get only data for a new field
00042 !
00043 ! !INPUT/OUTPUT PARAMETERS:
00044 !
00045       Type (Enddef_search), Intent (InOut) :: search
00046 
00047 !     Data on the points to be searched
00048 
00049 !
00050 ! !OUTPUT PARAMETERS:
00051 !
00052       integer, Intent (Out)               :: recv_req (ndim_3d+2, search%npart)
00053 
00054 !     Receive requests of non-blocking receives
00055 
00056       integer, Intent (Out)               :: ierror
00057 
00058 !     Returns the error code of PSMILe_Recv_req_coords_real;
00059 !             ierror = 0 : No error
00060 !             ierror > 0 : Severe error
00061 !
00062 ! !LOCAL VARIABLES
00063 !
00064 !     ... for the intersections
00065 !
00066       integer                      :: i, ipart
00067       integer                      :: dim (ndim_3d)
00068       integer                      :: len_mask
00069 !
00070 !     ... for error parameters
00071 !
00072       integer, parameter           :: nerrp = 3
00073       integer                      :: ierrp (nerrp)
00074 !
00075 ! !DESCRIPTION:
00076 !
00077 ! Subroutine "PSMILe_Recv_req_coords_real" receives the subgrid requested
00078 ! from the sending process. The subgrid is sent after the corresponding
00079 ! request was received by routine "PSMILe_Get_intersect".
00080 !
00081 !
00082 ! !REVISION HISTORY:
00083 !
00084 !   Date      Programmer   Description
00085 ! ----------  ----------   -----------
00086 ! 03.07.10    H. Ritzdorf  created
00087 !
00088 !EOP
00089 !----------------------------------------------------------------------
00090 !
00091 ! $Id: psmile_recv_req_coords_real.F90 2325 2010-04-21 15:00:07Z valcke $
00092 ! $Author: valcke $
00093 !
00094    Character(len=len_cvs_string), save :: mycvs = 
00095        '$Id: psmile_recv_req_coords_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00096 !
00097 !----------------------------------------------------------------------
00098 !
00099 !  Initialization
00100 !
00101       ierror = 0
00102 
00103 #ifdef VERBOSE
00104       print 9990, trim(ch_id), sender, search%npart, receive_mask
00105 
00106       call psmile_flushstd
00107 #endif /* VERBOSE */
00108 !
00109       if (new_search) then
00110          Allocate (search%search_real (ndim_3d, search%npart), stat = ierror)
00111          if ( ierror /= 0 ) then
00112             ierrp (1) = search%npart * ndim_3d
00113             ierror = PRISM_Error_Alloc
00114             call psmile_error ( ierror, 'search%search_real', &
00115                   ierrp, 1, __FILE__, __LINE__ )
00116             return
00117          endif
00118       endif
00119 !
00120 ! ====================================================================
00121 !     Get dimensions of arrays to be allocated and
00122 !     allocate the arrays.
00123 ! ====================================================================
00124 
00125       do ipart = 1, search%npart
00126 !
00127       select case ( search%grid_type )
00128 
00129 ! -----------------------------------------------------------------------
00130 !         Regular in all directions
00131 ! -----------------------------------------------------------------------
00132 !
00133       case (PRISM_Reglonlatvrt)
00134 !
00135          dim(:) = search%range (2,:, ipart) - search%range (1,:, ipart) + 1
00136 
00137          len_mask = dim (1) * dim (2) * dim (3)
00138 !
00139 ! -----------------------------------------------------------------------
00140 !       Gauss Reduced grid in lonlat direction
00141 !         Regular in vertical direction
00142 !         For now search%range(*,2,*) is 1!
00143 ! -----------------------------------------------------------------------
00144 !
00145       case (PRISM_Gaussreduced_regvrt)
00146 
00147          dim(1) = (search%range(2,1, ipart) - search%range(1,1, ipart) + 1) * &
00148                   (search%range(2,2, ipart) - search%range(1,2, ipart) + 1)
00149          dim(2) = dim(1)
00150          dim(3) =  search%range(2,3, ipart) - search%range(1,3, ipart) + 1
00151 
00152          len_mask = dim (1) * dim (3)
00153 !
00154 ! -----------------------------------------------------------------------
00155 !       Irregular in lonlat   direction
00156 !         Regular in vertical direction
00157 ! -----------------------------------------------------------------------
00158 !
00159       case (PRISM_Irrlonlat_regvrt)
00160 
00161          dim(1) = (search%range(2,1, ipart) - search%range(1,1, ipart) + 1) * &
00162                   (search%range(2,2, ipart) - search%range(1,2, ipart) + 1)
00163          dim(2) = dim(1)
00164          dim(3) =  search%range(2,3, ipart) - search%range(1,3, ipart) + 1
00165 
00166          len_mask = dim (1) * dim (3)
00167 
00168 ! -----------------------------------------------------------------------
00169 !       Irregular in lonlat   and vertical direction
00170 ! -----------------------------------------------------------------------
00171 !
00172       case (PRISM_Irrlonlatvrt)
00173 
00174          dim(1) = (search%range(2,1, ipart) - search%range(1,1, ipart) + 1) * &
00175                   (search%range(2,2, ipart) - search%range(1,2, ipart) + 1) * &
00176                   (search%range(2,3, ipart) - search%range(1,3, ipart) + 1)
00177          dim(2) = dim(1)
00178          dim(3) = dim(1)
00179 
00180          len_mask = dim (1)
00181 !
00182 ! -----------------------------------------------------------------------
00183 !     Error: unsupported grid type
00184 ! -----------------------------------------------------------------------
00185 !
00186       case DEFAULT
00187 !
00188           ierrp (1) = search%grid_type
00189 
00190           ierror = PRISM_Error_Internal
00191 
00192           call psmile_error ( ierror, 'unsupported grid generation type', &
00193                               ierrp, 1, __FILE__, __LINE__ )
00194 
00195       end select
00196 !
00197       search%dims(:, ipart) = dim (:)
00198 !
00199 !  ... Allocate arrays
00200 !
00201       if (new_search) then
00202          do i = 1, ndim_3d
00203 
00204          Allocate (search%search_real(i, ipart)%vector(dim(i)), STAT = ierror)
00205          if ( ierror > 0 ) then
00206             ierrp (1) = ierror
00207             ierrp (2) = dim(i)
00208             ierror = PRISM_Error_Alloc
00209             call psmile_error ( ierror, 'search%search_real(i)%vector', &
00210                                 ierrp, 2, __FILE__, __LINE__ )
00211             return
00212          endif
00213 
00214          end do
00215       endif
00216 !
00217 !===> ... Post receive
00218 !         !!! Das ist unschoen mit den vielen sends/receives.
00219 !         !!! Was schoeneres faellt mir aber im Moment nicht ein.
00220 !
00221 !     if (sender /= psmile_rank) then
00222 !     endif
00223 
00224          do i = 1, ndim_3d
00225 
00226          call MPI_Irecv (search%search_real(i, ipart)%vector, dim(i), &
00227                          MPI_REAL, &
00228                          sender, tag, comm_psmile, recv_req (i, ipart), ierror)
00229          if (ierror /= MPI_SUCCESS) then
00230             ierrp (1) = ierror
00231             ierrp (2) = sender
00232             ierrp (3) = tag
00233             ierror = PRISM_Error_Recv
00234 
00235             call psmile_error (ierror, 'MPI_Irecv', &
00236                                ierrp, 3, __FILE__, __LINE__ )
00237             return
00238          endif
00239 !
00240          end do
00241 !
00242 !===> Post receive for the mask if required
00243 !
00244       if (receive_mask) then
00245 
00246          if (new_search .or. .not. Associated (search%search_mask(ipart)%vector)) &
00247             then
00248 
00249             Allocate (search%search_mask(ipart)%vector(len_mask), STAT = ierror)
00250             if ( ierror > 0 ) then
00251                ierrp (1) = ierror
00252                ierrp (2) = len_mask
00253                ierror = PRISM_Error_Alloc
00254                call psmile_error ( ierror, 'search%search_mask(ipart)%vector', &
00255                                    ierrp, 2, __FILE__, __LINE__ )
00256                return
00257             endif
00258          endif
00259 
00260          call MPI_Irecv (search%search_mask(ipart)%vector, len_mask, &
00261                          MPI_LOGICAL, sender, tag, comm_psmile, &
00262                          recv_req (ndim_3d+1, ipart), ierror)
00263          if (ierror /= MPI_SUCCESS) then
00264             ierrp (1) = ierror
00265             ierrp (2) = sender
00266             ierrp (3) = tag
00267             ierror = PRISM_Error_Recv
00268 
00269             call psmile_error (ierror, 'MPI_Irecv of mask', &
00270                                ierrp, 3, __FILE__, __LINE__ )
00271             return
00272          endif
00273 
00274       else
00275 !
00276 !Rene: WARNING: We have to keep search_mask here when this call was preceded
00277 !      psmile_recv_req_corners_real which may have received a mask already.
00278 !      What shall we do if the cell-based search required the mask but
00279 !      the point-based search does not require a mask? We need more information
00280 !      from the target to make a decision, something like a no-mask flag. We could
00281 !      initiate a new search : new_search == .true. ???
00282 !
00283 !      Workaraound for the user in case:
00284 !        - only point-based search is required for a particular grid
00285 !        - some fields for this grid have a PRISM_Undefined mask
00286 !        ->  define a mask with all values flagged as .true. for those fields
00287 !        
00288 !
00289          if ( Associated(search%search_mask) ) then
00290             if (.not. new_search .and. Associated (search%search_mask(ipart)%vector)) &
00291                  then
00292 
00293 !rr               Deallocate (search%search_mask(ipart)%vector)
00294             endif
00295 
00296 !rr               Nullify    (search%search_mask(ipart)%vector)
00297 !rr               if ( ipart == search%npart ) Nullify (search%search_mask)
00298 
00299          endif
00300 
00301          recv_req (ndim_3d+1, ipart) = MPI_REQUEST_NULL
00302 
00303       endif
00304       end do ! ipart
00305 !
00306 !===> All done
00307 !
00308 #ifdef VERBOSE
00309       print 9980, trim(ch_id), ierror
00310 
00311       call psmile_flushstd
00312 #endif /* VERBOSE */
00313 !
00314 !  Formats:
00315 !
00316 9990 format (1x, a, ': PSMILe_Recv_req_coords_real: ', &
00317                     'sender =', i4, ', npart =', i3, ', receive_mask', l2)
00318 9980 format (1x, a, ': PSMILe_Recv_req_coords_real: eof ', &
00319                     'ierror =', i4)
00320 
00321       end subroutine PSMILe_Recv_req_coords_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1