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 3248 2011-06-23 13:03:19Z coquart $
00104 ! $Author: coquart $
00105 !
00106    Character(len=len_cvs_string), save :: mycvs = 
00107        '$Id: psmile_recv_req_subgrid.F90 3248 2011-06-23 13:03:19Z coquart $'
00108 !
00109 !----------------------------------------------------------------------
00110 !
00111 !  Initialization
00112 !
00113 #ifdef VERBOSE
00114       print 9990, trim(ch_id), msg_intersections%src_comp_id,              &
00115                                msg_intersections%field_info%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%field_info%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%search_data%range (2, ndim_3d, npart), &
00133                    search%search_data%shape (2, ndim_3d, npart), &
00134                    search%search_data%dim_size (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%search_data%range(:,:,n) = msg_intersections%intersections(n)%intersection
00157          end do
00158 !
00159          search%search_data%shape (1:2, 1:ndim_3d, 1:npart) = &
00160             search%search_data%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%search_data%npart       = npart
00168          search%method_type = method_type
00169          search%search_data%datatype    = datatype
00170 !
00171          sgrid = msg_intersections%intersections(1)%tgt_all_extents_grid_id
00172          search%search_data%grid_type   = &
00173             all_comp_infos(msg_intersections%all_comp_infos_comp_idx)%all_extent_infos(sgrid)%grid_type
00174 !
00175          search%len_msg     = nd_msgint
00176          search%sender      = sender
00177          if (Associated(search%boundary_cell)) Deallocate(search%boundary_cell)
00178 !
00179 ! HUHU: Temporary fix; some "shape_is_for_corners(ndim_3d)" variable should
00180 !       become part of search structure
00181 !
00182          shape_req_corners = (msg_intersections%field_info%requires_conserv_remap == PSMILe_conserv2D .or.   &
00183                               msg_intersections%field_info%requires_conserv_remap == PSMILe_conserv3D) .and. &
00184                              search%search_data%grid_type /= PRISM_Gridless
00185       endif
00186 !
00187       if (recv_mask) then
00188          if (new_search .or. .not. Associated (search%search_mask)) then
00189 
00190             Allocate (search%search_mask (npart), stat = ierror)
00191             if ( ierror /= 0 ) then
00192                ierrp (1) = npart
00193                ierror = PRISM_Error_Alloc
00194                call psmile_error ( ierror, 'search%search_mask', &
00195                      ierrp, 1, __FILE__, __LINE__ )
00196                return
00197             endif
00198             do n = 1, npart
00199                      Nullify (search%search_mask(n)%vector )
00200             enddo
00201          endif
00202 
00203       else
00204 
00205          if (.not. new_search .and. Associated (search%search_mask)) then
00206 !
00207 !Rene We have to keep search_mask here as this may have been transmitted
00208 !     already with the corners for conservative remapping'
00209 !
00210 !rr            Deallocate (search%search_mask)
00211          endif
00212 
00213 !rr         Nullify (search%search_mask)
00214          recv_req (ndim_3d+1, 1:npart) = MPI_REQUEST_NULL
00215 
00216       endif
00217 !
00218 !======================================================================
00219 !     Receive coordinates and mask for
00220 !======================================================================
00221 !
00222       if (search%search_data%grid_type == PRISM_Gridless) then
00223 !
00224 !----------------------------------------------------------------------
00225 !     ... a gridless grid (only the mask is sent)
00226 !
00227 !     TODO: PSMILe_Recv_req_mask in PSMILe_Recv_req_coords verwenden !
00228 !           Dazu muss aber die Reihenfolge geaendert werden.
00229 !           Deswegen wird auch noch der Move gemacht.
00230 !----------------------------------------------------------------------
00231 !
00232            call psmile_recv_req_mask (sender, tag, search,             &
00233                                       recv_req, recv_mask, new_search, &
00234                                       ierror)
00235 
00236 !          Move real request to the beginning
00237 
00238            temp = recv_req (1,1)
00239            recv_req (1,1) = recv_req (ndim_3d+1,1)
00240            recv_req (ndim_3d+1,1) = temp
00241 !
00242 !----------------------------------------------------------------------
00243 !     Conservative Remapping
00244 !     Search range is always the same as originally defined
00245 !----------------------------------------------------------------------
00246 !
00247       else if ( msg_intersections%field_info%requires_conserv_remap == PSMILe_conserv2D .or. &
00248                 msg_intersections%field_info%requires_conserv_remap == PSMILe_conserv3D ) then
00249 !
00250          if (datatype == MPI_REAL) then
00251 
00252 !        ... Real datatype
00253 
00254             call psmile_recv_req_corners_real (sender, tag, search, &
00255                                                recv_req, recv_mask, new_search, &
00256                                                ierror)
00257 
00258          else if (datatype == MPI_DOUBLE_PRECISION) then
00259 
00260 !        ... Double datatype
00261 
00262             call psmile_recv_req_corners_dble (sender, tag, search, &
00263                                                recv_req, recv_mask, new_search, &
00264                                                ierror)
00265 
00266 #if defined ( PRISM_QUAD_TYPE )
00267          else if (datatype == MPI_REAL16) then
00268 
00269 !        ... Quadruple  datatype
00270 
00271             call psmile_recv_req_corner_quad (sender, tag, search, &
00272                                               recv_req, recv_mask, new_search, &
00273                                               ierror)
00274 #endif
00275 
00276          else
00277 !
00278 !           Unknown data type
00279 !
00280             ierrp (1) = datatype
00281             ierror = PRISM_Error_Internal
00282             call psmile_error ( ierror, 'unsupported data type', &
00283                                 ierrp, 1, __FILE__, __LINE__ )
00284          endif
00285 
00286       else
00287 !
00288 !----------------------------------------------------------------------
00289 !     ... Point Method
00290 !         Recv range must be (temporarily) decreased by one
00291 !         if original range is defined for corners.
00292 !----------------------------------------------------------------------
00293 !
00294          if ( shape_req_corners .and. &
00295               search%search_data%grid_type /= PRISM_Gaussreduced_regvrt) then
00296 !
00297 !  Decrease dimensions by one as we have to switch from corners to points
00298 !  HUHU: This must depend on PSMILe_conserv2D/PSMILe_conserv3D
00299 !
00300             search%search_data%range (2, 1:ndim_2d, :) = search%search_data%range (2, 1:ndim_2d, :) - 1
00301          endif
00302 !
00303          select case ( search%method_type )
00304 
00305          case (PSMILe_PointMethod)
00306 !
00307             if (datatype == MPI_REAL) then
00308 
00309 !           ... Real datatype
00310 
00311                call psmile_recv_req_coords_real (sender, tag, search, &
00312                                                  recv_req, recv_mask, new_search, &
00313                                                  ierror)
00314 
00315             else if (datatype == MPI_DOUBLE_PRECISION) then
00316 
00317 !           ... Double datatype
00318 
00319                call psmile_recv_req_coords_dble (sender, tag, search, &
00320                                                  recv_req, recv_mask, new_search, &
00321                                                  ierror)
00322 
00323 #if defined ( PRISM_QUAD_TYPE )
00324             else if (datatype == MPI_REAL16) then
00325 
00326 !           ... Quadruple datatype
00327 
00328                call psmile_recv_req_coords_quad (sender, tag, search, &
00329                                                  recv_req, recv_mask, new_search, &
00330                                                  ierror)
00331 
00332 #endif
00333 
00334             else
00335 !
00336 !           ... Unknown data type
00337 !
00338                ierrp (1) = datatype
00339                ierror = PRISM_Error_Internal
00340                call psmile_error ( ierror, 'unsupported data type', &
00341                                    ierrp, 1, __FILE__, __LINE__ )
00342             endif
00343 
00344 !rr We do not need to increase the range as we have to decrease in psmle_search_donor_irreg2 again.
00345 
00346             if ( shape_req_corners .and. &
00347                  search%search_data%grid_type /= PRISM_Gaussreduced_regvrt ) then
00348 
00349                search%search_data%range (2, 1:ndim_2d, :) = search%search_data%range (2, 1:ndim_2d, :) + 1
00350             endif
00351 !
00352 !----------------------------------------------------------------------
00353 !        ... Vector method
00354 !----------------------------------------------------------------------
00355 !
00356          case (PSMILe_VectorPointMethod)
00357 
00358             ierrp (1) = method_type
00359             ierror = PRISM_Error_Internal
00360             call psmile_error ( ierror, 'Vector Method is currently not supported', &
00361                                 ierrp, 1, __FILE__, __LINE__ )
00362 !
00363 !----------------------------------------------------------------------
00364 !        ... Subgrid method
00365 !----------------------------------------------------------------------
00366 !
00367          case (PSMILe_SubgridMethod)
00368 
00369             ierrp (1) = method_type
00370             ierror = PRISM_Error_Internal
00371             call psmile_error ( ierror, 'Subgrid Method is currently not supported', &
00372                                 ierrp, 1, __FILE__, __LINE__ )
00373 !
00374 !----------------------------------------------------------------------
00375 !     ... Unknown method
00376 !----------------------------------------------------------------------
00377 !
00378          case default
00379 
00380             ierrp (1) = method_type
00381             ierror = PRISM_Error_Internal
00382             call psmile_error ( ierror, 'unsupported method type', &
00383                                 ierrp, 1, __FILE__, __LINE__ )
00384 
00385          end select ! case ( search%method_type )
00386 !
00387       endif
00388 !
00389 !===> All done
00390 !
00391 #ifdef VERBOSE
00392       print 9980, trim(ch_id), msg_intersections%src_comp_id, sender, ierror
00393 
00394       call psmile_flushstd
00395 #endif /* VERBOSE */
00396 !
00397 !  Formats:
00398 !
00399 9990 format (1x, a, ': psmile_recv_req_subgrid: comp_id =', i3, &
00400                     '; method id =', i4, '; sender =', i4, '; tag =', i4)
00401 9980 format (1x, a, ': psmile_recv_req_subgrid: comp_id =', i3, &
00402                     '; eof sender =', i3, ', ierror =', i4)
00403 
00404       end subroutine PSMILe_Recv_req_subgrid

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1