psmile_recv_req_subgrid.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_subgrid
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_recv_req_subgrid (msg_intersections, &
00012                                           sender, tag, search, recv_req, &
00013                                           new_search, ierror)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM_constants
00018 !
00019       use PSMILe, dummy_interface => PSMILe_Recv_req_subgrid
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       Type (enddef_msg_intersections), Intent (In) :: msg_intersections
00026 
00027 !     The request for a subgrid.
00028 !     (cf. description in include file "psmile_common.F90"
00029 
00030       Integer, Intent (In)                :: sender
00031 
00032 !     Specifies the sender of the grid.
00033 
00034       Integer, Intent (In)                :: tag
00035 
00036 !     Specifies the message tag used
00037 
00038       Logical, Intent (In)                :: new_search
00039 
00040 !     Option:
00041 !     new_search = .true.  : New search; fill search structure
00042 !     new_search = .false. : Get only data for a new field
00043 !
00044 ! !INPUT/OUTPUT PARAMETERS:
00045 !
00046       Type (Enddef_search), Intent (InOut) :: search
00047 
00048 !     Data on the points to be searched
00049 !
00050       integer, Intent (InOut)               :: recv_req (ndim_3d+2, *)
00051 
00052 !     Receive requests of non-blocking receives!
00053 !
00054 ! !OUTPUT PARAMETERS:
00055 !
00056       integer, Intent (Out)               :: ierror
00057 
00058 !     Returns the error code of PSMILe_Recv_req_subgrid;
00059 !             ierror = 0 : No error
00060 !             ierror > 0 : Severe error
00061 !
00062 ! !LOCAL VARIABLES
00063 !
00064       integer                      :: sgrid
00065 !
00066 !     ... number of intersections
00067 !
00068       integer                      :: n, npart
00069 !
00070 !     ... for the point information
00071 !
00072       logical, save                :: shape_req_corners
00073       integer                      :: method_type
00074       integer                      :: datatype
00075 !
00076       logical                      :: recv_mask
00077 !
00078 !     ... for receive requests
00079 !
00080       integer                      :: temp
00081 !
00082 !     ... for error parameters
00083 !
00084       integer, parameter           :: nerrp = 1
00085       integer                      :: ierrp (nerrp)
00086 !
00087 ! !DESCRIPTION:
00088 !
00089 ! Subroutine "PSMILe_Recv_req_subgrid" receives the subgrid requested from
00090 ! the sending process. The subgrid was sent after the corresponding
00091 ! request was received by routine "PSMILe_Get_intersect".
00092 !
00093 !
00094 ! !REVISION HISTORY:
00095 !
00096 !   Date      Programmer   Description
00097 ! ----------  ----------   -----------
00098 ! 03.07.10    H. Ritzdorf  created
00099 !
00100 !EOP
00101 !----------------------------------------------------------------------
00102 !
00103 ! $Id: psmile_recv_req_subgrid.F90 2781 2010-11-26 15:07:14Z hanke $
00104 ! $Author: hanke $
00105 !
00106    Character(len=len_cvs_string), save :: mycvs = 
00107        '$Id: psmile_recv_req_subgrid.F90 2781 2010-11-26 15:07:14Z hanke $'
00108 !
00109 !----------------------------------------------------------------------
00110 !
00111 !  Initialization
00112 !
00113 #ifdef VERBOSE
00114       print 9990, trim(ch_id), msg_intersections%src_comp_id,         &
00115                                msg_intersections%first_tgt_method_id, &
00116                                sender, tag
00117       print *, 'new_search = ', new_search
00118       print *, 'ndim_3d = ', ndim_3d 
00119       call psmile_flushstd
00120 #endif /* VERBOSE */
00121 !
00122       method_type = msg_intersections%method_type
00123       datatype    = msg_intersections%method_datatype
00124       recv_mask   = msg_intersections%tgt_mask_id /= PRISM_UNDEFINED
00125 !
00126       npart       = msg_intersections%num_parts
00127 !
00128       if (new_search) then
00129 !
00130 !===> ... Allocate vectors
00131 !
00132          Allocate (search%range (2, ndim_3d, npart), &
00133                    search%shape (2, ndim_3d, npart), &
00134                    search%dims  (   ndim_3d, npart), &
00135                    search%msgint(nd_msgint), stat = ierror)
00136 
00137          if ( ierror /= 0 ) then
00138             ierrp (1) = 5 * ndim_3d * npart + nd_msgint
00139             ierror = PRISM_Error_Alloc
00140             call psmile_error ( ierror, 'search%{range,shape,dims,msgint}', &
00141                   ierrp, 1, __FILE__, __LINE__ )
00142             return
00143          endif
00144 
00145          Nullify (search%search_mask)
00146 !
00147 !===> ... Copy range
00148 !         Note: The initial range is always sufficent to store 
00149 !               data for corner search which requires a larger dimension
00150 !               than point search.
00151 !               HUHU: Should be probably reallocated in order
00152 !                     to speed up the memory access ?!?
00153 !                     See also other HUHU remarks.
00154 !
00155          do n = 1, npart
00156             search%range(:,:,n) = msg_intersections%intersections(n)%intersection
00157          end do
00158 !
00159          search%shape (1:2, 1:ndim_3d, 1:npart) = &
00160             search%range (1:2, 1:ndim_3d, 1:npart)
00161 !
00162          call psmile_copy_msg_intersections (search%msg_intersections, msg_intersections)
00163          call psmile_pack_msg_intersections (search%msg_intersections, search%msgint)
00164 !
00165 !===> ... Copy further values
00166 !
00167          search%npart       = npart
00168          search%method_type = method_type
00169          search%datatype    = datatype
00170 !
00171          sgrid = msg_intersections%intersections(1)%tgt_all_extents_grid_id
00172          search%grid_type   = all_comp_infos(msg_intersections%all_comp_infos_comp_idx)%all_extent_infos(2, sgrid)
00173 !
00174          search%len_msg     = nd_msgint
00175          search%sender      = sender
00176          if (Associated(search%boundary_cell)) Deallocate(search%boundary_cell)
00177 !
00178 ! HUHU: Temporary fix; some "shape_is_for_corners(ndim_3d)" variable should
00179 !       become part of search structure
00180 !
00181          shape_req_corners = (msg_intersections%requires_conserv_remap == PSMILe_conserv2D .or.   &
00182                               msg_intersections%requires_conserv_remap == PSMILe_conserv3D) .and. &
00183                              search%grid_type /= PRISM_Gridless
00184       endif
00185 !
00186       if (recv_mask) then
00187          if (new_search .or. .not. Associated (search%search_mask)) then
00188 
00189             Allocate (search%search_mask (npart), stat = ierror)
00190             if ( ierror /= 0 ) then
00191                ierrp (1) = npart
00192                ierror = PRISM_Error_Alloc
00193                call psmile_error ( ierror, 'search%search_mask', &
00194                      ierrp, 1, __FILE__, __LINE__ )
00195                return
00196             endif
00197             do n = 1, npart
00198                      Nullify (search%search_mask(n)%vector )
00199             enddo
00200          endif
00201 
00202       else
00203 
00204          if (.not. new_search .and. Associated (search%search_mask)) then
00205 !
00206 !Rene We have to keep search_mask here as this may have been transmitted
00207 !     already with the corners for conservative remapping'
00208 !
00209 !rr            Deallocate (search%search_mask)
00210          endif
00211 
00212 !rr         Nullify (search%search_mask)
00213          recv_req (ndim_3d+1, 1:npart) = MPI_REQUEST_NULL
00214 
00215       endif
00216 !
00217 !======================================================================
00218 !     Receive coordinates and mask for
00219 !======================================================================
00220 !
00221       if (search%grid_type == PRISM_Gridless) then
00222 !
00223 !----------------------------------------------------------------------
00224 !     ... a gridless grid (only the mask is sent)
00225 !
00226 !     TODO: PSMILe_Recv_req_mask in PSMILe_Recv_req_coords verwenden !
00227 !           Dazu muss aber die Reihenfolge geaendert werden.
00228 !           Deswegen wird auch noch der Move gemacht.
00229 !----------------------------------------------------------------------
00230 !
00231            call psmile_recv_req_mask (sender, tag, search,             &
00232                                       recv_req, recv_mask, new_search, &
00233                                       ierror)
00234 
00235 !          Move real request to the beginning
00236 
00237            temp = recv_req (1,1)
00238            recv_req (1,1) = recv_req (ndim_3d+1,1)
00239            recv_req (ndim_3d+1,1) = temp
00240 !
00241 !----------------------------------------------------------------------
00242 !     Conservative Remapping
00243 !     Search range is always the same as originally defined
00244 !----------------------------------------------------------------------
00245 !
00246       else if ( msg_intersections%requires_conserv_remap == PSMILe_conserv2D .or. &
00247                 msg_intersections%requires_conserv_remap == PSMILe_conserv3D ) then
00248 !
00249          if (datatype == MPI_REAL) then
00250 
00251 !        ... Real datatype
00252 
00253             call psmile_recv_req_corners_real (sender, tag, search, &
00254                                                recv_req, recv_mask, new_search, &
00255                                                ierror)
00256 
00257          else if (datatype == MPI_DOUBLE_PRECISION) then
00258 
00259 !        ... Double datatype
00260 
00261             call psmile_recv_req_corners_dble (sender, tag, search, &
00262                                                recv_req, recv_mask, new_search, &
00263                                                ierror)
00264 
00265 #if defined ( PRISM_QUAD_TYPE )
00266          else if (datatype == MPI_REAL16) then
00267 
00268 !        ... Quadruple  datatype
00269 
00270             call psmile_recv_req_corner_quad (sender, tag, search, &
00271                                               recv_req, recv_mask, new_search, &
00272                                               ierror)
00273 #endif
00274 
00275          else
00276 !
00277 !           Unknown data type
00278 !
00279             ierrp (1) = datatype
00280             ierror = PRISM_Error_Internal
00281             call psmile_error ( ierror, 'unsupported data type', &
00282                                 ierrp, 1, __FILE__, __LINE__ )
00283          endif
00284 
00285       else
00286 !
00287 !----------------------------------------------------------------------
00288 !     ... Point Method
00289 !         Recv range must be (temporarily) decreased by one
00290 !         if original range is defined for corners.
00291 !----------------------------------------------------------------------
00292 !
00293          if ( shape_req_corners .and. &
00294               search%grid_type /= PRISM_Gaussreduced_regvrt) then
00295 !
00296 !  Decrease dimensions by one as we have to switch from corners to points
00297 !  HUHU: This must depend on PSMILe_conserv2D/PSMILe_conserv3D
00298 !
00299             search%range (2, 1:ndim_2d, :) = search%range (2, 1:ndim_2d, :) - 1
00300          endif
00301 !
00302          select case ( search%method_type )
00303 
00304          case (PSMILe_PointMethod)
00305 !
00306             if (datatype == MPI_REAL) then
00307 
00308 !           ... Real datatype
00309 
00310                call psmile_recv_req_coords_real (sender, tag, search, &
00311                                                  recv_req, recv_mask, new_search, &
00312                                                  ierror)
00313 
00314             else if (datatype == MPI_DOUBLE_PRECISION) then
00315 
00316 !           ... Double datatype
00317 
00318                call psmile_recv_req_coords_dble (sender, tag, search, &
00319                                                  recv_req, recv_mask, new_search, &
00320                                                  ierror)
00321 
00322 #if defined ( PRISM_QUAD_TYPE )
00323             else if (datatype == MPI_REAL16) then
00324 
00325 !           ... Quadruple datatype
00326 
00327                call psmile_recv_req_coords_quad (sender, tag, search, &
00328                                                  recv_req, recv_mask, new_search, &
00329                                                  ierror)
00330 
00331 #endif
00332 
00333             else
00334 !
00335 !           ... Unknown data type
00336 !
00337                ierrp (1) = datatype
00338                ierror = PRISM_Error_Internal
00339                call psmile_error ( ierror, 'unsupported data type', &
00340                                    ierrp, 1, __FILE__, __LINE__ )
00341             endif
00342 
00343 !rr We do not need to increase the range as we have to decrease in psmle_search_donor_irreg2 again.
00344 
00345             if ( shape_req_corners .and. &
00346                  search%grid_type /= PRISM_Gaussreduced_regvrt ) then
00347 
00348                search%range (2, 1:ndim_2d, :) = search%range (2, 1:ndim_2d, :) + 1
00349             endif
00350 !
00351 !----------------------------------------------------------------------
00352 !        ... Vector method
00353 !----------------------------------------------------------------------
00354 !
00355          case (PSMILe_VectorPointMethod)
00356 
00357             ierrp (1) = method_type
00358             ierror = PRISM_Error_Internal
00359             call psmile_error ( ierror, 'Vector Method is currently not supported', &
00360                                 ierrp, 1, __FILE__, __LINE__ )
00361 !
00362 !----------------------------------------------------------------------
00363 !        ... Subgrid method
00364 !----------------------------------------------------------------------
00365 !
00366          case (PSMILe_SubgridMethod)
00367 
00368             ierrp (1) = method_type
00369             ierror = PRISM_Error_Internal
00370             call psmile_error ( ierror, 'Subgrid Method is currently not supported', &
00371                                 ierrp, 1, __FILE__, __LINE__ )
00372 !
00373 !----------------------------------------------------------------------
00374 !     ... Unknown method
00375 !----------------------------------------------------------------------
00376 !
00377          case default
00378 
00379             ierrp (1) = method_type
00380             ierror = PRISM_Error_Internal
00381             call psmile_error ( ierror, 'unsupported method type', &
00382                                 ierrp, 1, __FILE__, __LINE__ )
00383 
00384          end select ! case ( search%method_type )
00385 !
00386       endif
00387 !
00388 !===> All done
00389 !
00390 #ifdef VERBOSE
00391       print 9980, trim(ch_id), msg_intersections%src_comp_id, sender, ierror
00392 
00393       call psmile_flushstd
00394 #endif /* VERBOSE */
00395 !
00396 !  Formats:
00397 !
00398 9990 format (1x, a, ': psmile_recv_req_subgrid: comp_id =', i3, &
00399                     '; method id =', i4, '; sender =', i4, '; tag =', i4)
00400 9980 format (1x, a, ': psmile_recv_req_subgrid: comp_id =', i3, &
00401                     '; eof sender =', i3, ', ierror =', i4)
00402 
00403       end subroutine PSMILe_Recv_req_subgrid

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1