psmile_get_info_index.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, NEC Europe Ltd., London, UK.
00003 ! Copyright 2010, DKRZ, Hamburg, Germany.
00004 ! All rights reserved. Use is subject to OASIS4 license terms.
00005 !-----------------------------------------------------------------------
00006 !BOP
00007 !
00008 ! !ROUTINE: PSMILe_Get_info_index
00009 !
00010 ! !INTERFACE:
00011 
00012       subroutine psmile_get_info_index ( method_id, request, index, ierror )
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_Get_info_index
00019 !
00020       use psmile_reallocate, only : psmile_realloc
00021 !
00022       implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026       Integer, Intent (in) :: method_id 
00027 
00028 !     Handle to method for which an index should be returned.
00029 
00030       Integer, Intent (in) :: request
00031 
00032 !     Type of index which should be returned
00033 !     request = 1 (send_coupler_index) : Get index for send to coupler
00034 !     request = 2 (send_direct_index)  : Get index for direct send
00035 !     request = 3 (recv_coupler_index) : Get index for recv from coupler
00036 !     request = 4 (recv_direct_index)  : Get index for direct recv
00037 !     request = 5 (send_appl_index)    : Get index for send to application
00038 !
00039 ! !OUTPUT PARAMETERS:
00040 !
00041       integer, Intent (Out) :: index
00042 
00043 !     Index corresponding to request
00044 
00045       integer, Intent (Out) :: ierror
00046 
00047 !     Returns the error code of PSMILe_Get_info_index;
00048 !             ierror = 0 : No error
00049 !             ierror > 0 : Severe error
00050 !
00051 !
00052 ! !LOCAL PARAMETERS:
00053 !
00054       Integer, Parameter :: send_coupler_index = 1
00055       Integer, Parameter :: send_direct_index  = 2
00056       Integer, Parameter :: recv_coupler_index = 3
00057       Integer, Parameter :: recv_direct_index  = 4
00058       Integer, Parameter :: send_appl_index    = 5
00059 !
00060       Integer, Parameter :: new_alloc = 8
00061 !
00062 ! !LOCAL VARIABLES
00063 !
00064       Type (Method), Pointer :: mp
00065 !
00066       Integer, Save :: unique_appl_msg_id = 1
00067 !
00068 !
00069 ! !DESCRIPTION:
00070 !
00071 ! Subroutine "PSMILe_Get_info_index" returns an index for an 
00072 ! send or receive info within a method with method id "method_id".
00073 !
00074 ! !REVISION HISTORY:
00075 !
00076 !   Date      Programmer    Description
00077 ! ----------  -----------   -----------
00078 ! 01.12.03    H. Ritzdorf   created
00079 ! 02.12.10    M. Hanke      Introduced usage of psmile_realloc
00080 !
00081 !EOP
00082 !----------------------------------------------------------------------
00083 !
00084 ! $Id: psmile_get_info_index.F90 2799 2010-12-03 14:33:35Z hanke $
00085 ! $Author: hanke $
00086 !
00087    Character(len=len_cvs_string), save :: mycvs = 
00088        '$Id: psmile_get_info_index.F90 2799 2010-12-03 14:33:35Z hanke $'
00089 !
00090 !----------------------------------------------------------------------
00091 
00092 #ifdef VERBOSE
00093       print 9990, trim(ch_id), method_id
00094 
00095       call psmile_flushstd
00096 #endif /* VERBOSE */
00097 
00098 !   Initialization
00099 !
00100       ierror = 0
00101 !
00102 #ifdef PRISM_ASSERTION
00103       if (method_id < 1 .or. &
00104           method_id > Number_of_Methods_allocated ) then
00105          call psmile_assert ( __FILE__, __LINE__, "illegal method id" )
00106          return
00107       endif
00108 #endif
00109 !
00110       mp => Methods (method_id)
00111 !
00112 !     Get index depending on the request
00113 !
00114       select case (request)
00115 
00116 !-----------------------------------------------------------------------
00117 !     Get index for send to coupler process
00118 !-----------------------------------------------------------------------
00119 
00120       case (send_coupler_index)
00121 
00122          if (mp%n_send_info_coupler == 0) nullify (mp%send_infos_coupler)
00123 
00124          mp%n_send_info_coupler = mp%n_send_info_coupler + 1
00125          index = mp%n_send_info_coupler
00126 
00127          if (mp%n_send_info_coupler > mp%n_alloc_send_coupler) then
00128 
00129            mp%n_alloc_send_coupler = mp%n_alloc_send_coupler + new_alloc
00130            mp%send_infos_coupler => psmile_realloc (mp%send_infos_coupler, &
00131                                                     mp%n_alloc_send_coupler)
00132          endif
00133 
00134 !-----------------------------------------------------------------------
00135 !     Get index for direct send to MPI process of another application
00136 !-----------------------------------------------------------------------
00137 
00138       case (send_direct_index)
00139 
00140          if (mp%n_send_info_direct == 0) nullify (mp%send_infos_direct)
00141 
00142          mp%n_send_info_direct = mp%n_send_info_direct + 1
00143          index = mp%n_send_info_direct
00144 
00145          if (mp%n_send_info_direct > mp%n_alloc_send_direct) then
00146 
00147             mp%n_alloc_send_direct = mp%n_alloc_send_direct + new_alloc
00148             mp%send_infos_direct => psmile_realloc (mp%send_infos_direct, &
00149                                                     mp%n_alloc_send_direct)
00150          endif
00151 
00152 !-----------------------------------------------------------------------
00153 !     Get index for receive from a coupler process
00154 !-----------------------------------------------------------------------
00155 !
00156       case (recv_coupler_index)
00157 
00158          if (mp%n_recv_info_coupler == 0) nullify (mp%recv_infos_coupler)
00159 
00160          mp%n_recv_info_coupler = mp%n_recv_info_coupler + 1
00161          index = mp%n_recv_info_coupler
00162 
00163          if (mp%n_recv_info_coupler > mp%n_alloc_recv_coupler) then
00164 
00165            mp%n_alloc_recv_coupler = mp%n_alloc_recv_coupler + new_alloc
00166            mp%recv_infos_coupler => psmile_realloc (mp%recv_infos_coupler, &
00167                                                     mp%n_alloc_recv_coupler)
00168          endif
00169 
00170 !-----------------------------------------------------------------------
00171 !     Get index for direct receive from MPI process of another application
00172 !-----------------------------------------------------------------------
00173 
00174       case (recv_direct_index)
00175 
00176          if (mp%n_recv_info_direct == 0) nullify (mp%recv_infos_direct)
00177 
00178          mp%n_recv_info_direct = mp%n_recv_info_direct + 1
00179          index = mp%n_recv_info_direct
00180 
00181          if (mp%n_recv_info_direct > mp%n_alloc_recv_direct) then
00182 
00183            mp%n_alloc_recv_direct = mp%n_alloc_recv_direct + new_alloc
00184            mp%recv_infos_direct => psmile_realloc (mp%recv_infos_direct, &
00185                                                    mp%n_alloc_recv_direct)
00186          endif
00187 
00188 !-----------------------------------------------------------------------
00189 !     Get index for send to a MPI process of same application
00190 !     (used to send data of global search)
00191 !-----------------------------------------------------------------------
00192 
00193       case (send_appl_index)
00194 
00195          if (mp%n_send_info_appl == 0) nullify (mp%send_infos_appl)
00196 
00197          mp%n_send_info_appl = mp%n_send_info_appl + 1
00198          index = mp%n_send_info_appl
00199 
00200          if (mp%n_send_info_appl > mp%n_alloc_send_appl) then
00201 
00202            mp%n_alloc_send_appl = mp%n_send_info_appl + new_alloc
00203            mp%send_infos_appl => psmile_realloc (mp%send_infos_appl, &
00204                                                  mp%n_alloc_send_appl)
00205          endif
00206 
00207          mp%send_infos_appl(index)%msg_id = unique_appl_msg_id
00208          unique_appl_msg_id = unique_appl_msg_id + 1
00209 !
00210 !-----------------------------------------------------------------------
00211 !     Invalid request
00212 !-----------------------------------------------------------------------
00213 
00214       case DEFAULT
00215 
00216 !        ... Error: invalid request
00217 
00218          ierror = PRISM_Error_Internal
00219 
00220          call psmile_error ( ierror, 'Invalid request', &
00221                              (/method_id, request/), 2, &
00222                              __FILE__, __LINE__ )
00223          return
00224 
00225       end select
00226 
00227 #ifdef VERBOSE
00228       print 9980, trim(ch_id), ierror, index
00229 
00230       call psmile_flushstd
00231 #endif /* VERBOSE */
00232 !
00233 !  Formats:
00234 !     
00235 
00236 #ifdef VERBOSE
00237 
00238 9990 format (1x, a, ': psmile_get_info_index: method_id', i3)
00239 9980 format (1x, a, ': psmile_get_info_index: eof ierror =', i3, &
00240                     '; index =', i4)
00241 
00242 #endif /* VERBOSE */
00243 !
00244       end subroutine PSMILe_Get_info_index

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1