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, 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 
00026       Type (enddef_msg_extra), Intent (In)  :: msg_extra
00027 
00028 !
00029 !     Derived type containing the request for an extra search.
00030 !
00031       Integer,            Intent (In)       :: sender
00032 !
00033 !     Rank of sending process in communicator "comm_psmile".
00034 !
00035 ! !OUTPUT PARAMETERS:
00036 !
00037       Integer,            Intent (Out)      :: ierror
00038 
00039 !     Returns the error code of PSMILe_Get_locations_3d;
00040 !             ierror = 0 : No error
00041 !             ierror > 0 : Severe error
00042 !
00043 ! !LOCAL PARAMETERS
00044 !
00045 !  TOL must become input parameter
00046 !
00047       Double Precision, parameter           :: tol = 1d-6
00048 !
00049 ! !LOCAL VARIABLES
00050 
00051       Integer                               :: n_buf, n_send, msg_extra_buf(msg_extra_size)
00052 #ifdef PRISM_ASSERTION
00053       Integer                               :: n_elem
00054 #endif
00055 !
00056       Type (Enddef_global_search)           :: search_global
00057 
00058 !     ... for point-to-point communication
00059 
00060       Integer                               :: status (MPI_STATUS_SIZE)
00061 !
00062 !     ... for error handling
00063 !
00064       Integer, Parameter                    :: nerrp = 2
00065       Integer                               :: ierrp (nerrp)
00066 !
00067 ! !DESCRIPTION:
00068 !
00069 ! Subroutine "PSMILe_Enddef_action_extra" performs the actions in order to
00070 ! search extra points for interpolation. The request was sent by process
00071 ! "sender".
00072 ! The data is sent by routine "PSMILe_global_search_..." in the process
00073 ! which has initually searched the data.
00074 !
00075 ! !REVISION HISTORY:
00076 !
00077 !   Date      Programmer   Description
00078 ! ----------  ----------   -----------
00079 ! 21.07.05    H. Ritzdorf  created
00080 !
00081 !EOP
00082 !----------------------------------------------------------------------
00083 !
00084 !  $Id: psmile_enddef_action_extra.F90 3112 2011-04-07 15:03:18Z hanke $
00085 !  $Author: hanke $
00086 !
00087    Character(len=len_cvs_string), save :: mycvs = 
00088        '$Id: psmile_enddef_action_extra.F90 3112 2011-04-07 15:03:18Z hanke $'
00089 !
00090 !----------------------------------------------------------------------
00091 !
00092 !  Initialization
00093 !
00094 #ifdef VERBOSE
00095       print 9990, trim(ch_id), sender
00096 
00097       call psmile_flushstd
00098 #endif /* VERBOSE */
00099 !
00100       n_send = msg_extra%len_int_data
00101       n_buf  = msg_extra%len_coord_data
00102 !
00103 #ifdef PRISM_ASSERTION
00104       if (n_buf <= 0) then
00105          print *, 'n_buf =', n_buf
00106          call psmile_assert (__FILE__, __LINE__, &
00107               "n_buf should be > 0")
00108       endif
00109 
00110       if (n_send <= 0) then
00111          print *, 'n_send =', n_send
00112          call psmile_assert (__FILE__, __LINE__, &
00113               "n_send should be > 0")
00114       endif
00115 #endif
00116 !
00117 !======================================================================
00118 !  Allocate integer buffer in order to receive the integer data
00119 !======================================================================
00120 !
00121       Allocate (search_global%ibuf (n_send), STAT = ierror)
00122 
00123       if ( ierror > 0 ) then
00124          ierrp (1) = ierror
00125          ierrp (2) = n_send
00126 
00127          ierror = PRISM_Error_Alloc
00128          call psmile_error ( ierror, 'search_global%ibuf', &
00129                              ierrp, 2, __FILE__, __LINE__ )
00130          return
00131       endif
00132 !
00133       call MPI_Recv (search_global%ibuf, n_send, MPI_INTEGER, &
00134                      sender, exttag, comm_psmile, status, ierror)
00135       if ( ierror /= MPI_SUCCESS ) then
00136          ierrp (1) = ierror
00137          ierror = PRISM_Error_MPI
00138 
00139          call psmile_error ( ierror, 'MPI_Recv(search_global%ibuf)', &
00140                              ierrp, 1, __FILE__, __LINE__ )
00141          return
00142       endif
00143 !
00144 #ifdef PRISM_ASSERTION
00145       call MPI_Get_elements (status, MPI_INTEGER, n_elem, ierror)
00146       if (ierror == MPI_SUCCESS) then
00147          if (n_elem /= n_send) then
00148             ierrp (1) = n_elem
00149             ierrp (2) = n_send
00150             ierror = PRISM_Error_Internal
00151 
00152             call psmile_error ( ierror, 'Incorrect message length: found/expected', &
00153                                 ierrp, 2, __FILE__, __LINE__ )
00154             return
00155 
00156          endif
00157       endif
00158 #endif
00159 !
00160 !======================================================================
00161 !     Receive real data on colntrol volumes and coordinates
00162 !     to be searched.
00163 !     Switch depending on the datatype to be searched extra
00164 !======================================================================
00165 !
00166       select case (msg_extra%datatype)
00167 
00168       case (PRISM_Real)
00169 !
00170 !       Allocate real buffer
00171 !       in order to receive the real data
00172 !
00173          Allocate (search_global%rbuf (n_buf), STAT = ierror)
00174 !
00175          if ( ierror > 0 ) then
00176             ierrp (1) = ierror
00177             ierrp (2) = n_buf
00178 
00179             ierror = PRISM_Error_Alloc
00180             call psmile_error ( ierror, 'search_global%rbuf', &
00181                                 ierrp, 2, __FILE__, __LINE__ )
00182             return
00183          endif
00184 !
00185          call MPI_Recv (search_global%rbuf, n_buf, MPI_REAL, &
00186                         sender, exttag, comm_psmile, status, ierror)
00187          if ( ierror /= MPI_SUCCESS ) then
00188             ierrp (1) = ierror
00189             ierror = PRISM_Error_MPI
00190 
00191             call psmile_error ( ierror, 'MPI_Recv(search_global%rbuf)', &
00192                                 ierrp, 1, __FILE__, __LINE__ )
00193             return
00194          endif
00195 !
00196 #ifdef PRISM_ASSERTION
00197          call MPI_Get_elements (status, MPI_REAL, n_elem, ierror)
00198          if (ierror == MPI_SUCCESS) then
00199             if (n_elem /= n_buf) then
00200                ierrp (1) = n_elem
00201                ierrp (2) = n_buf
00202                ierror = PRISM_Error_Internal
00203 
00204                call psmile_error ( ierror, 'Incorrect message length: found/expected', &
00205                                    ierrp, 2, __FILE__, __LINE__ )
00206                return
00207 
00208             endif
00209          endif
00210 #endif
00211 
00212       case (PRISM_Double_Precision)
00213 !
00214 !       Allocate double precision buffer
00215 !       in order to receive the double precision data
00216 !
00217          Allocate (search_global%dbuf (n_buf), STAT = ierror)
00218 !
00219          if ( ierror > 0 ) then
00220             ierrp (1) = ierror
00221             ierrp (2) = n_buf
00222 
00223             ierror = PRISM_Error_Alloc
00224             call psmile_error ( ierror, 'search_global%dbuf', &
00225                                 ierrp, 2, __FILE__, __LINE__ )
00226             return
00227          endif
00228 !
00229          call MPI_Recv (search_global%dbuf, n_buf, MPI_DOUBLE_PRECISION, &
00230                         sender, exttag, comm_psmile, status, ierror)
00231          if ( ierror /= MPI_SUCCESS ) then
00232             ierrp (1) = ierror
00233             ierror = PRISM_Error_MPI
00234 
00235             call psmile_error ( ierror, 'MPI_Recv(search_global%dbuf)', &
00236                                 ierrp, 1, __FILE__, __LINE__ )
00237             return
00238          endif
00239 !
00240 #ifdef PRISM_ASSERTION
00241          call MPI_Get_elements (status, MPI_DOUBLE_PRECISION, n_elem, ierror)
00242          if (ierror == MPI_SUCCESS) then
00243             if (n_elem /= n_buf) then
00244                ierrp (1) = n_elem
00245                ierrp (2) = n_buf
00246                ierror = PRISM_Error_Internal
00247 
00248                call psmile_error ( ierror, 'Incorrect message length: found/expected', &
00249                                    ierrp, 2, __FILE__, __LINE__ )
00250                return
00251 
00252             endif
00253          endif
00254 #endif
00255 
00256 #if defined ( PRISM_QUAD_TYPE )
00257       case (TODO_QUAD)
00258 !
00259 !       Allocate quadruple precision buffer
00260 !       in order to receive the quadruple precision data
00261 !
00262          Allocate (search_global%qbuf (n_buf), STAT = ierror)
00263 !
00264          if ( ierror > 0 ) then
00265             ierrp (1) = ierror
00266             ierrp (2) = n_buf
00267 
00268             ierror = PRISM_Error_Alloc
00269             call psmile_error ( ierror, 'search_global%qbuf', &
00270                                 ierrp, 2, __FILE__, __LINE__ )
00271             return
00272          endif
00273 
00274         TODO: This has to be implemented
00275 #endif
00276 
00277       case default
00278          ierror = PRISM_Error_Internal
00279          ierrp (1) = msg_extra%datatype
00280 
00281          call psmile_error ( ierror, 'Unknown datatype', &
00282                              ierrp, 1, __FILE__, __LINE__ )
00283          return
00284 
00285       end select
00286 !
00287 !======================================================================
00288 !     Search Data in local 
00289 !======================================================================
00290 !
00291       search_global%sender = sender
00292       search_global%msg_extra = msg_extra
00293 
00294       call psmile_search_donor_extra (search_global, tol, ierror)
00295       if (ierror > 0) return
00296 
00297 #ifdef FOO
00298 !======================================================================
00299 !     Send results back to requesting process
00300 !     Sending process has already posted a corresonding receive
00301 !======================================================================
00302 !
00303       call psmile_pack_msg_extra(msg_extra, msg_extra_buf)
00304 
00305       call MPI_Send (msg_extra_buf, msg_extra_size, MPI_INTEGER, &
00306                      sender, rexttag, comm_psmile, ierror)
00307 
00308       if (ierror /= MPI_SUCCESS) then
00309          ierrp (1) = ierror
00310          ierror = PRISM_Error_MPI
00311 
00312          call psmile_error ( ierror, 'MPI_Send(msg_extra_buf)', &
00313                              ierrp, 1, __FILE__, __LINE__ )
00314          return
00315       endif
00316 #endif
00317 !
00318 !===> All done
00319 !
00320 #ifdef VERBOSE
00321       print 9980, trim(ch_id), ierror
00322 
00323       call psmile_flushstd
00324 #endif /* VERBOSE */
00325 !
00326 !  Formats:
00327 !
00328 #ifdef VERBOSE
00329 
00330 9990 format (1x, a, ': psmile_enddef_action_extra: sender ', i6)
00331 9980 format (1x, a, ': psmile_enddef_action_extra: eof ierror =', i3)
00332 
00333 #endif /* VERBOSE */
00334 
00335       end subroutine PSMILe_Enddef_action_extra

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1