psmile_enddef_action_extra.F90

Go to the documentation of this file.
00001 !
00002 !-----------------------------------------------------------------------
00003 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00004 ! All rights reserved. Use is subject to OASIS4 license terms.
00005 !-----------------------------------------------------------------------
00006 !BOP
00007 !
00008 ! !ROUTINE: PSMILe_Enddef_action_extra
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_enddef_action_extra (msg_extra, nd_msg, sender, &
00013                         ierror)
00014 !
00015 ! !USES:
00016 !
00017       use PRISM
00018 !
00019       use PSMILe, dummy_interface => PSMILe_Enddef_action_extra
00020 
00021       implicit none
00022 !
00023 ! !INPUT PARAMETERS:
00024 !
00025       Integer,            Intent (In)  :: nd_msg
00026 
00027 !     Number of data items in vector "msg_extra".
00028 
00029       Integer,            Intent (In)  :: msg_extra (nd_msg)
00030 
00031 !     Integer vector containing the request for an extra search.
00032 !
00033 !  extra_msg (1)  = Type of request
00034 !                   PSMILE_Trilinear : Search additional points for
00035 !                                      trilinear interpolation.
00036 !  extra_msg (2)  = Datatype
00037 !  extra_msg (3)  = Length of buffer "ibuf" containing the integer data.
00038 !  extra_msg (4)  = Length of buffer "buf" containing the
00039 !                   real/double precision data.
00040 !  extra_msg (5)  = Global comp id of component
00041 !  extra_msg (6)  = Global var id
00042 !  extra_msg (7)  = Grid Type of sending grid
00043 !                   (PRISM_Reglonlatvrt, ...)
00044 !  extra_msg (8)  = Number of real/double data items per control
00045 !                   volume sent
00046 !  extra_msg (9)  = Is data on global offset vailable ?
00047 !  extra_msg (10) = I-Index of global offset (partition)
00048 !  extra_msg (11) = J-Index of global offset (partition)
00049 !  extra_msg (12) = K-Index of global offset (partition)
00050 !
00051       Integer,            Intent (In)  :: sender
00052 !
00053 !     Rank of sending process in communicator "comm_psmile".
00054 !
00055 ! !OUTPUT PARAMETERS:
00056 !
00057       Integer,            Intent (Out) :: ierror
00058 
00059 !     Returns the error code of PSMILe_Get_locations_3d;
00060 !             ierror = 0 : No error
00061 !             ierror > 0 : Severe error
00062 !
00063 ! !LOCAL PARAMETERS
00064 !
00065 !  TOL must become input parameter
00066 !
00067       Double Precision, parameter     :: tol = 1d-6
00068 !
00069 ! !LOCAL VARIABLES
00070 
00071       Integer                         :: n_buf, n_send
00072 #ifdef PRISM_ASSERTION
00073       Integer                         :: n_elem
00074 #endif
00075 !
00076       Type (Enddef_global_search)     :: search_global
00077 
00078 !     ... for point-to-point communication
00079 
00080       Integer                         :: status (MPI_STATUS_SIZE)
00081 !
00082 !     ... for error handling
00083 !
00084       Integer, Parameter              :: nerrp = 2
00085       Integer                         :: ierrp (nerrp)
00086 !
00087 ! !DESCRIPTION:
00088 !
00089 ! Subroutine "PSMILe_Enddef_action_extra" performs the actions in order to
00090 ! search extra points for interpolation. The request was sent by process
00091 ! "sender".
00092 ! The data is sent by routine "PSMILe_global_search_..." in the process
00093 ! which has initually searched the data.
00094 !
00095 ! !REVISION HISTORY:
00096 !
00097 !   Date      Programmer   Description
00098 ! ----------  ----------   -----------
00099 ! 21.07.05    H. Ritzdorf  created
00100 !
00101 !EOP
00102 !----------------------------------------------------------------------
00103 !
00104 !  $Id: psmile_enddef_action_extra.F90 2082 2009-10-21 13:31:19Z hanke $
00105 !  $Author: hanke $
00106 !
00107    Character(len=len_cvs_string), save :: mycvs = 
00108        '$Id: psmile_enddef_action_extra.F90 2082 2009-10-21 13:31:19Z hanke $'
00109 !
00110 !----------------------------------------------------------------------
00111 !
00112 !  Initialization
00113 !
00114 #ifdef VERBOSE
00115       print 9990, trim(ch_id), sender
00116 
00117       call psmile_flushstd
00118 #endif /* VERBOSE */
00119 !
00120       n_send = msg_extra (3)
00121       n_buf  = msg_extra (4)
00122 !
00123 #ifdef PRISM_ASSERTION
00124       if (n_buf <= 0) then
00125          print *, 'n_buf =', n_buf
00126          call psmile_assert (__FILE__, __LINE__, &
00127               "n_buf should be > 0")
00128       endif
00129 
00130       if (n_send <= 0) then
00131          print *, 'n_send =', n_send
00132          call psmile_assert (__FILE__, __LINE__, &
00133               "n_send should be > 0")
00134       endif
00135 #endif
00136 !
00137 !======================================================================
00138 !  Allocate integer buffer in order to receive the integer data
00139 !======================================================================
00140 !
00141       Allocate (search_global%ibuf (n_send), STAT = ierror)
00142 
00143       if ( ierror > 0 ) then
00144          ierrp (1) = ierror
00145          ierrp (2) = n_send
00146 
00147          ierror = PRISM_Error_Alloc
00148          call psmile_error ( ierror, 'search_global%ibuf', &
00149                              ierrp, 2, __FILE__, __LINE__ )
00150          return
00151       endif
00152 !
00153       call MPI_Recv (search_global%ibuf, n_send, MPI_INTEGER, &
00154                      sender, exttag, comm_psmile, status, ierror)
00155       if ( ierror /= MPI_SUCCESS ) then
00156          ierrp (1) = ierror
00157          ierror = PRISM_Error_MPI
00158 
00159          call psmile_error ( ierror, 'MPI_Recv(search_global%ibuf)', &
00160                              ierrp, 1, __FILE__, __LINE__ )
00161          return
00162       endif
00163 !
00164 #ifdef PRISM_ASSERTION
00165       call MPI_Get_elements (status, MPI_INTEGER, n_elem, ierror)
00166       if (ierror == MPI_SUCCESS) then
00167          if (n_elem /= n_send) then
00168             ierrp (1) = n_elem
00169             ierrp (2) = n_send
00170             ierror = PRISM_Error_Internal
00171 
00172             call psmile_error ( ierror, 'Incorrect message length: found/expected', &
00173                                 ierrp, 2, __FILE__, __LINE__ )
00174             return
00175 
00176          endif
00177       endif
00178 #endif
00179 !
00180 !======================================================================
00181 !     Receive real data on colntrol volumes and coordinates
00182 !     to be searched.
00183 !     Switch depending on the datatype to be searched extra
00184 !======================================================================
00185 !
00186       select case (msg_extra(2))
00187 
00188       case (PRISM_Real)
00189 !
00190 !       Allocate real buffer
00191 !       in order to receive the real data
00192 !
00193          Allocate (search_global%rbuf (n_buf), STAT = ierror)
00194 !
00195          if ( ierror > 0 ) then
00196             ierrp (1) = ierror
00197             ierrp (2) = n_buf
00198 
00199             ierror = PRISM_Error_Alloc
00200             call psmile_error ( ierror, 'search_global%rbuf', &
00201                                 ierrp, 2, __FILE__, __LINE__ )
00202             return
00203          endif
00204 !
00205          call MPI_Recv (search_global%rbuf, n_buf, MPI_REAL, &
00206                         sender, exttag, comm_psmile, status, ierror)
00207          if ( ierror /= MPI_SUCCESS ) then
00208             ierrp (1) = ierror
00209             ierror = PRISM_Error_MPI
00210 
00211             call psmile_error ( ierror, 'MPI_Recv(search_global%rbuf)', &
00212                                 ierrp, 1, __FILE__, __LINE__ )
00213             return
00214          endif
00215 !
00216 #ifdef PRISM_ASSERTION
00217          call MPI_Get_elements (status, MPI_REAL, n_elem, ierror)
00218          if (ierror == MPI_SUCCESS) then
00219             if (n_elem /= n_buf) then
00220                ierrp (1) = n_elem
00221                ierrp (2) = n_buf
00222                ierror = PRISM_Error_Internal
00223 
00224                call psmile_error ( ierror, 'Incorrect message length: found/expected', &
00225                                    ierrp, 2, __FILE__, __LINE__ )
00226                return
00227 
00228             endif
00229          endif
00230 #endif
00231 
00232       case (PRISM_Double_Precision)
00233 !
00234 !       Allocate double precision buffer
00235 !       in order to receive the double precision data
00236 !
00237          Allocate (search_global%dbuf (n_buf), STAT = ierror)
00238 !
00239          if ( ierror > 0 ) then
00240             ierrp (1) = ierror
00241             ierrp (2) = n_buf
00242 
00243             ierror = PRISM_Error_Alloc
00244             call psmile_error ( ierror, 'search_global%dbuf', &
00245                                 ierrp, 2, __FILE__, __LINE__ )
00246             return
00247          endif
00248 !
00249          call MPI_Recv (search_global%dbuf, n_buf, MPI_DOUBLE_PRECISION, &
00250                         sender, exttag, comm_psmile, status, ierror)
00251          if ( ierror /= MPI_SUCCESS ) then
00252             ierrp (1) = ierror
00253             ierror = PRISM_Error_MPI
00254 
00255             call psmile_error ( ierror, 'MPI_Recv(search_global%dbuf)', &
00256                                 ierrp, 1, __FILE__, __LINE__ )
00257             return
00258          endif
00259 !
00260 #ifdef PRISM_ASSERTION
00261          call MPI_Get_elements (status, MPI_DOUBLE_PRECISION, n_elem, ierror)
00262          if (ierror == MPI_SUCCESS) then
00263             if (n_elem /= n_buf) then
00264                ierrp (1) = n_elem
00265                ierrp (2) = n_buf
00266                ierror = PRISM_Error_Internal
00267 
00268                call psmile_error ( ierror, 'Incorrect message length: found/expected', &
00269                                    ierrp, 2, __FILE__, __LINE__ )
00270                return
00271 
00272             endif
00273          endif
00274 #endif
00275 
00276 #if defined ( PRISM_QUAD_TYPE )
00277       case (TODO_QUAD)
00278 !
00279 !       Allocate quadruple precision buffer
00280 !       in order to receive the quadruple precision data
00281 !
00282          Allocate (search_global%qbuf (n_buf), STAT = ierror)
00283 !
00284          if ( ierror > 0 ) then
00285             ierrp (1) = ierror
00286             ierrp (2) = n_buf
00287 
00288             ierror = PRISM_Error_Alloc
00289             call psmile_error ( ierror, 'search_global%qbuf', &
00290                                 ierrp, 2, __FILE__, __LINE__ )
00291             return
00292          endif
00293 
00294         TODO: This has to be implemented
00295 #endif
00296 
00297       case default
00298          ierror = PRISM_Error_Internal
00299          ierrp (1) = msg_extra(2)
00300 
00301          call psmile_error ( ierror, 'Unknown datatype', &
00302                              ierrp, 1, __FILE__, __LINE__ )
00303          return
00304 
00305       end select
00306 !
00307 !======================================================================
00308 !     Search Data in local 
00309 !======================================================================
00310 !
00311       Allocate (search_global%msg_extra(1:nd_msg), STAT = ierror)
00312 !
00313       if ( ierror > 0 ) then
00314          ierrp (1) = ierror
00315          ierrp (2) = nd_msg
00316 
00317          ierror = PRISM_Error_Alloc
00318          call psmile_error ( ierror, 'search_global%msgint', &
00319                              ierrp, 2, __FILE__, __LINE__ )
00320          return
00321       endif
00322 !
00323       search_global%sender = sender
00324       search_global%msg_extra = msg_extra (1:nd_msg)
00325 !     search_global%len_msg = nd_msg
00326 !
00327       call psmile_search_donor_extra (search_global, tol, ierror)
00328       if (ierror > 0) return
00329 
00330 #ifdef FOO
00331 !======================================================================
00332 !     Send results back to requesting process
00333 !     Sending process has already posted a corresonding receive
00334 !======================================================================
00335 !
00336       call MPI_Send (msg_extra, nd_msg, MPI_INTEGER, &
00337                      sender, rexttag, comm_psmile, ierror)
00338 
00339       if (ierror /= MPI_SUCCESS) then
00340          ierrp (1) = ierror
00341          ierror = PRISM_Error_MPI
00342 
00343          call psmile_error ( ierror, 'MPI_Send(msg_extra)', &
00344                              ierrp, 1, __FILE__, __LINE__ )
00345          return
00346       endif
00347 #endif
00348 !
00349 !===> All done
00350 !
00351 #ifdef VERBOSE
00352       print 9980, trim(ch_id), ierror
00353 
00354       call psmile_flushstd
00355 #endif /* VERBOSE */
00356 !
00357 !  Formats:
00358 !
00359 #ifdef VERBOSE
00360 
00361 9990 format (1x, a, ': psmile_enddef_action_extra: sender ', i6)
00362 9980 format (1x, a, ': psmile_enddef_action_extra: eof ierror =', i3)
00363 
00364 #endif /* VERBOSE */
00365 
00366       end subroutine PSMILe_Enddef_action_extra

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1