psmile_get_exch_index.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_Get_exch_index
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_get_exch_index ( var_id, task_id, request, index, ierror )
00012 !
00013 ! !USES:
00014 !
00015       use PRISM_constants
00016 !
00017       use PSMILe, dummy_interface => PSMILe_Get_exch_index
00018 !
00019       implicit none
00020 !
00021 ! !INPUT PARAMETERS:
00022 !
00023       Integer, Intent (in) :: var_id
00024 
00025 !     Handle to method for which an index should be returned.
00026 
00027       Integer, Intent (in) :: task_id
00028 
00029 !     Task Id relevant for the send task
00030 
00031       Integer, Intent (in) :: request
00032 
00033 !     Type of index which should be returned
00034 !     request = 1 (send_coupler_index) : Get index for send to coupler
00035 !     request = 2 (send_direct_index)  : Get index for direct send
00036 !     request = 3 (recv_coupler_index) : Get index for recv from coupler
00037 !     request = 4 (recv_direct_index)  : Get index for direct recv
00038 !     request = 5 (send_appl_index)    : Get index for send to application
00039 !
00040 ! !OUTPUT PARAMETERS:
00041 !
00042       integer, Intent (Out) :: index
00043 
00044 !     Index corresponding to request
00045 
00046       integer, Intent (Out) :: ierror
00047 
00048 !     Returns the error code of PSMILe_Get_exch_index;
00049 !             ierror = 0 : No error
00050 !             ierror > 0 : Severe error
00051 !
00052 !
00053 ! !LOCAL PARAMETERS:
00054 !
00055       Integer, Parameter :: send_coupler_index = 1
00056       Integer, Parameter :: send_direct_index  = 2
00057       Integer, Parameter :: recv_coupler_index = 3
00058       Integer, Parameter :: recv_direct_index  = 4
00059       Integer, Parameter :: send_appl_index    = 5
00060 !
00061       Integer, Parameter :: new_alloc = 2
00062 !
00063 ! !LOCAL VARIABLES
00064 !
00065       Type (Taskout_type), Pointer :: fieldout
00066       Type (Taskin_type ), Pointer :: fieldin
00067 
00068       integer             :: new_dim
00069 
00070       Type (Send_field_information), Pointer :: sendinfo (:)
00071       Type (Recv_field_information), Pointer :: recvinfo (:)
00072 
00073 !     ... for error handling
00074 
00075       Integer, parameter  :: nerrp = 2
00076       Integer             :: ierrp (nerrp)
00077 !
00078 !
00079 ! !DESCRIPTION:
00080 !
00081 ! Subroutine "PSMILe_Get_exch_index" returns an index for an 
00082 ! send or receive info within a field with field id "var_id".
00083 !
00084 ! !REVISION HISTORY:
00085 !
00086 !   Date      Programmer    Description
00087 ! ----------  -----------   -----------
00088 ! 01.12.03    H. Ritzdorf   created
00089 !
00090 !EOP
00091 !----------------------------------------------------------------------
00092 !
00093 ! $Id: psmile_get_exch_index.F90 2325 2010-04-21 15:00:07Z valcke $
00094 ! $Author: valcke $
00095 !
00096    Character(len=len_cvs_string), save :: mycvs = 
00097     '$Id: psmile_get_exch_index.F90 2325 2010-04-21 15:00:07Z valcke $'
00098 !
00099 !----------------------------------------------------------------------
00100 
00101 #ifdef VERBOSE
00102       print 9990, trim(ch_id), var_id
00103 
00104       call psmile_flushstd
00105 #endif /* VERBOSE */
00106 
00107 !   Initialization
00108 !
00109       ierror = 0
00110 !
00111 #ifdef PRISM_ASSERTION
00112       if (var_id < 1 .or. &
00113           var_id > Number_of_Fields_allocated ) then
00114          print *, 'var id', var_id
00115          call psmile_assert ( __FILE__, __LINE__, "illegal var id" )
00116          return
00117       endif
00118 #endif
00119 !
00120       if ( request == send_coupler_index .or. &
00121            request == send_appl_index    .or. &
00122            request == send_direct_index ) then
00123          fieldout => Fields(var_id)%Taskout(task_id)
00124 !
00125       else if ( request == recv_coupler_index .or. &
00126                 request == recv_direct_index ) then
00127          fieldin  => Fields(var_id)%Taskin
00128       endif
00129 !
00130 !     Get index depending on the request
00131 !
00132       select case (request)
00133 
00134 !-----------------------------------------------------------------------
00135 !     Get index for send to coupler process
00136 !-----------------------------------------------------------------------
00137 
00138       case (send_coupler_index)
00139 
00140         fieldout%n_send_coupler = fieldout%n_send_coupler + 1
00141         index = fieldout%n_send_coupler
00142 
00143         if (fieldout%n_send_coupler > fieldout%n_alloc_send_coupler) then
00144            new_dim = fieldout%n_alloc_send_coupler + new_alloc
00145 !
00146            Allocate (sendinfo (new_dim), STAT = ierror)
00147            if (ierror > 0) then
00148               ierrp (1) = ierror
00149               ierrp (2) = new_dim
00150               ierror = PRISM_Error_Alloc
00151 
00152               call psmile_error ( ierror, 'sendinfo', &
00153                                   ierrp, 2, __FILE__, __LINE__ )
00154               return
00155            endif
00156 !
00157            if (fieldout%n_alloc_send_coupler > 0) then
00158               sendinfo (1:fieldout%n_alloc_send_coupler) = &
00159                  fieldout%send_coupler (1:fieldout%n_alloc_send_coupler)
00160 !
00161 !          De-allocate old vector
00162 !
00163               Deallocate (fieldout%send_coupler, STAT = ierror)
00164               if (ierror > 0) then
00165                  ierrp (1) = ierror
00166                  ierror = PRISM_Error_Dealloc
00167 
00168                  call psmile_error ( ierror, 'fieldout%send_coupler', &
00169                                      ierrp, 1, __FILE__, __LINE__ )
00170                  return
00171               endif
00172            endif
00173 !
00174            fieldout%send_coupler => sendinfo
00175            fieldout%n_alloc_send_coupler = new_dim
00176          endif
00177 
00178 !-----------------------------------------------------------------------
00179 !     Get index for direct send to MPI process of another application
00180 !-----------------------------------------------------------------------
00181 
00182       case (send_direct_index)
00183 
00184         fieldout%n_send_direct = fieldout%n_send_direct + 1
00185         index = fieldout%n_send_direct
00186 
00187         if (fieldout%n_send_direct > fieldout%n_alloc_send_direct) then
00188            new_dim = fieldout%n_alloc_send_direct + new_alloc
00189 !
00190            Allocate (sendinfo (new_dim), STAT = ierror)
00191            if (ierror > 0) then
00192               ierrp (1) = ierror
00193               ierrp (2) = new_dim
00194               ierror = PRISM_Error_Alloc
00195 
00196               call psmile_error ( ierror, 'sendinfo', &
00197                                   ierrp, 2, __FILE__, __LINE__ )
00198               return
00199            endif
00200 !
00201            if (fieldout%n_alloc_send_direct > 0) then
00202               sendinfo (1:fieldout%n_alloc_send_direct) = &
00203                  fieldout%send_direct (1:fieldout%n_alloc_send_direct)
00204 !
00205 !          De-allocate old vector
00206 !
00207               Deallocate (fieldout%send_direct, STAT = ierror)
00208               if (ierror > 0) then
00209                  ierrp (1) = ierror
00210                  ierror = PRISM_Error_Dealloc
00211 
00212                  call psmile_error ( ierror, 'fieldout%send_direct', &
00213                                      ierrp, 1, __FILE__, __LINE__ )
00214                  return
00215               endif
00216            endif
00217 !
00218            fieldout%send_direct => sendinfo
00219            fieldout%n_alloc_send_direct = new_dim
00220          endif
00221 !
00222 !-----------------------------------------------------------------------
00223 !     Get index for receive from a coupler process
00224 !-----------------------------------------------------------------------
00225 !
00226       case (recv_coupler_index)
00227 
00228         fieldin%n_recv_coupler = fieldin%n_recv_coupler + 1
00229         index = fieldin%n_recv_coupler
00230 
00231         if (fieldin%n_recv_coupler > fieldin%n_alloc_recv_coupler) then
00232            new_dim = fieldin%n_alloc_recv_coupler + new_alloc
00233 !
00234            Allocate (recvinfo (new_dim), STAT = ierror)
00235            if (ierror > 0) then
00236               ierrp (1) = ierror
00237               ierrp (2) = new_dim
00238               ierror = PRISM_Error_Alloc
00239 
00240               call psmile_error ( ierror, 'recvinfo', &
00241                                   ierrp, 2, __FILE__, __LINE__ )
00242               return
00243            endif
00244 !
00245            if (fieldin%n_alloc_recv_coupler > 0) then
00246               recvinfo (1:fieldin%n_alloc_recv_coupler) = &
00247                  fieldin%recv_coupler (1:fieldin%n_alloc_recv_coupler)
00248 !
00249 !          De-allocate old vector
00250 !
00251               Deallocate (fieldin%recv_coupler, STAT = ierror)
00252               if (ierror > 0) then
00253                  ierrp (1) = ierror
00254                  ierror = PRISM_Error_Dealloc
00255 
00256                  call psmile_error ( ierror, 'fieldin%recv_coupler', &
00257                                      ierrp, 1, __FILE__, __LINE__ )
00258                  return
00259               endif
00260            endif
00261 !
00262            fieldin%recv_coupler => recvinfo
00263            fieldin%n_alloc_recv_coupler = new_dim
00264          endif
00265 
00266 !-----------------------------------------------------------------------
00267 !     Get index for direct receive from MPI process of another application
00268 !-----------------------------------------------------------------------
00269 
00270       case (recv_direct_index)
00271 
00272         fieldin%n_recv_direct = fieldin%n_recv_direct + 1
00273         index = fieldin%n_recv_direct
00274 
00275         if (fieldin%n_recv_direct > fieldin%n_alloc_recv_direct) then
00276            new_dim = fieldin%n_alloc_recv_direct + new_alloc
00277 !
00278            Allocate (recvinfo (new_dim), STAT = ierror)
00279            if (ierror > 0) then
00280               ierrp (1) = ierror
00281               ierrp (2) = new_dim
00282               ierror = PRISM_Error_Alloc
00283 
00284               call psmile_error ( ierror, 'recvinfo', &
00285                                   ierrp, 2, __FILE__, __LINE__ )
00286               return
00287            endif
00288 !
00289            if (fieldin%n_alloc_recv_direct > 0) then
00290               recvinfo (1:fieldin%n_alloc_recv_direct) = &
00291                  fieldin%recv_direct (1:fieldin%n_alloc_recv_direct)
00292 !
00293 !          De-allocate old vector
00294 !
00295               Deallocate (fieldin%recv_direct, STAT = ierror)
00296               if (ierror > 0) then
00297                  ierrp (1) = ierror
00298                  ierror = PRISM_Error_Dealloc
00299 
00300                  call psmile_error ( ierror, 'fieldin%recv_direct', &
00301                                      ierrp, 1, __FILE__, __LINE__ )
00302                  return
00303               endif
00304            endif
00305 !
00306            fieldin%recv_direct => recvinfo
00307            fieldin%n_alloc_recv_direct = new_dim
00308          endif
00309 
00310 !-----------------------------------------------------------------------
00311 !     Get index for send to a MPI process of same application
00312 !     (used to send data of global search)
00313 !-----------------------------------------------------------------------
00314 
00315       case (send_appl_index)
00316 
00317         fieldout%n_send_appl = fieldout%n_send_appl + 1
00318         index = fieldout%n_send_appl
00319 
00320         if (fieldout%n_send_appl > fieldout%n_alloc_send_appl) then
00321            new_dim = fieldout%n_alloc_send_appl + new_alloc
00322 !
00323            Allocate (sendinfo (new_dim), STAT = ierror)
00324            if (ierror > 0) then
00325               ierrp (1) = ierror
00326               ierrp (2) = new_dim
00327               ierror = PRISM_Error_Alloc
00328 
00329               call psmile_error ( ierror, 'sendinfo', &
00330                                   ierrp, 2, __FILE__, __LINE__ )
00331               return
00332            endif
00333 !
00334            if (fieldout%n_alloc_send_appl > 0) then
00335               sendinfo (1:fieldout%n_alloc_send_appl) = &
00336                  fieldout%send_appl (1:fieldout%n_alloc_send_appl)
00337 !
00338 !          De-allocate old vector
00339 !
00340               Deallocate (fieldout%send_appl, STAT = ierror)
00341               if (ierror > 0) then
00342                  ierrp (1) = ierror
00343                  ierror = PRISM_Error_Dealloc
00344 
00345                  call psmile_error ( ierror, 'fieldout%send_appl', &
00346                                      ierrp, 1, __FILE__, __LINE__ )
00347                  return
00348               endif
00349            endif
00350 !
00351            fieldout%send_appl => sendinfo
00352            fieldout%n_alloc_send_appl = new_dim
00353          endif
00354 !
00355 !-----------------------------------------------------------------------
00356 !     Invalid request
00357 !-----------------------------------------------------------------------
00358 
00359       case DEFAULT
00360 
00361 !        ... Error: invalid request
00362 
00363          ierrp (1) = var_id
00364          ierrp (2) = request
00365 
00366          ierror = PRISM_Error_Internal
00367 
00368          call psmile_error ( ierror, 'Invalid request', &
00369                              ierrp, 2, __FILE__, __LINE__ )
00370          return
00371 
00372       end select
00373 
00374 #ifdef VERBOSE
00375       print 9980, trim(ch_id), ierror, index
00376 
00377       call psmile_flushstd
00378 #endif /* VERBOSE */
00379 !
00380 !  Formats:
00381 !     
00382 
00383 #ifdef VERBOSE
00384 
00385 9990 format (1x, a, ': psmile_get_exch_index: var_id', i3)
00386 9980 format (1x, a, ': psmile_get_exch_index: eof ierror =', i3, &
00387                     '; index =', i4)
00388 
00389 #endif /* VERBOSE */
00390 !
00391       end subroutine PSMILe_Get_exch_index

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1