psmile_get_face_ind_3d.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2008-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_face_ind_3d
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_get_face_ind_3d (search, extra_search,          &
00012                            send_info, len_cpl,                          &
00013                            send_mask, nreq, srcloc_ind, n_send,         &
00014                            ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM
00019 !
00020       use PSMILe, dummy_interface => PSMILe_Get_face_ind_3d
00021 
00022       Implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026       Type (Enddef_search),     Intent (In) :: search
00027 
00028 !     Info's on coordinates to be searched
00029 !
00030       Type (Extra_search_info), Intent (In) :: extra_search
00031 !
00032 !     Info's on locations which 
00033 !     (*) require global search or
00034 !     (*) extra neighbourhood search
00035 !
00036       Type (Send_information),  Intent (In) :: send_info
00037 !
00038 !     Send info's for sends to the coupler.
00039 !
00040       Integer,                  Intent (In) :: len_cpl (search%npart)
00041 
00042 !     Number of points to be sent to the coupler for each partition.
00043 
00044       Integer,                  Intent (In) :: nreq
00045 !
00046 !     Total number of extra points to be searched.
00047 !     Dimension of "send_mask".
00048 !
00049       Logical,                  Intent (In) :: send_mask (nreq)
00050 !
00051 !     Mask for the points (of "indices_req") to be sent.
00052 !
00053       Integer,                  Intent (In) :: n_send
00054 !
00055 !     Number of extra points to be sent and 
00056 !     number of indices to be returned.
00057 !
00058 ! !OUTPUT PARAMETERS:
00059 !
00060       Integer,                  Intent (Out) :: srcloc_ind (ndim_3d, n_send)
00061 !
00062 !     3d-indices of the locations to be sent.
00063 !
00064       Integer,                  Intent (Out) :: ierror
00065 
00066 !     Returns the error code of PSMILe_Get_face_ind_3d;
00067 !             ierror = 0 : No error
00068 !             ierror > 0 : Severe error
00069 !
00070 ! !LOCAL VARIABLES
00071 !
00072 !     ... loop variables
00073 !
00074       Integer                         :: i, j, n
00075 !
00076 !     ... for partitions
00077 !
00078       Integer                         :: ipart, nextra_prev
00079 !
00080       Integer, Pointer                :: indices_req (:)
00081       Integer, Pointer                :: len_req (:)
00082       Integer, Pointer                :: srcloc (:)
00083 !
00084 !     ... for error handling
00085 !
00086 !     Integer, Parameter              :: nerrp = 3
00087 !     Integer                         :: ierrp (nerrp)
00088 !
00089 ! !DESCRIPTION:
00090 !
00091 ! Subroutine "PSMILe_Get_face_ind_3d" computes the 3d indices of
00092 ! the locations to be sent to a "coinciding" (neighboured)
00093 ! process of the same component.
00094 !
00095 !
00096 ! !REVISION HISTORY:
00097 !
00098 !   Date      Programmer   Description
00099 ! ----------  ----------   -----------
00100 ! 02.02.05    H. Ritzdorf  created
00101 !
00102 !EOP
00103 !----------------------------------------------------------------------
00104 !
00105 !  $Id: psmile_get_face_ind_3d.F90 2933 2011-01-31 16:42:27Z hanke $
00106 !  $Author: hanke $
00107 !
00108    Character(len=len_cvs_string), save :: mycvs = 
00109        '$Id: psmile_get_face_ind_3d.F90 2933 2011-01-31 16:42:27Z hanke $'
00110 !
00111 !----------------------------------------------------------------------
00112 !
00113 !  Initialization
00114 !
00115 #ifdef VERBOSE
00116       print 9990, trim(ch_id)
00117 
00118       call psmile_flushstd
00119 #endif /* VERBOSE */
00120 !
00121       ierror = 0
00122 !
00123       len_req => extra_search%len_req
00124 !
00125 #ifdef PRISM_ASSERTION
00126       if (SUM(len_req(:)) > nreq) then
00127          print *, 'nreq, sum(len_req)', nreq, SUM(len_req(:))
00128          call psmile_assert ( __FILE__, __LINE__, &
00129                  'nreq is too small (inconsistent data)')
00130       endif
00131 !
00132       if (send_info%nvec /= 1) then
00133          print *, 'nvec', send_info%nvec
00134          call psmile_assert ( __FILE__, __LINE__, &
00135                  'Routine is designed for irregular case (3d vector)')
00136       endif
00137 #endif
00138 !
00139 !  For all parts "ipart"
00140 !
00141 !  indices_req = Indices of the extra points in entire list of all points
00142 !                to be searched (i.e. for all parts "ipart")
00143 !
00144       n = 0
00145       nextra_prev = 0
00146 !
00147 ! Note: Only 1 vector is generated for full 3d indices
00148 !       independent on number of partitions "search%npart".
00149 !       Thus, "nprev" is always 0 and independent on "len_cpl".
00150 !
00151       srcloc => send_info%srclocs(1,1)%vector
00152 !
00153       do ipart = 1, search%npart
00154 !
00155          if (len_req (ipart) > 0) then
00156             indices_req => extra_search%indices_req(ipart)%vector
00157 !
00158                do j = 1, len_req (ipart)
00159                if (send_mask(nextra_prev+j)) then
00160 !
00161 !  ... Compute index in srcloc
00162 !      "srcloc" is generated as array of dimension (ndim_3d, *).
00163 !      "indices_req(j) - nprev" is the index in the second dimension.
00164 !
00165                   n = n + 1
00166 !
00167                   i = (indices_req(j) - 1) * ndim_3d
00168 !
00169                   srcloc_ind (1, n) = srcloc (i+1)
00170                   srcloc_ind (2, n) = srcloc (i+2)
00171                   srcloc_ind (3, n) = srcloc (i+3)
00172                endif
00173                end do ! j
00174 !
00175             nextra_prev = nextra_prev + len_req (ipart)
00176          endif
00177 !
00178       end do ! ipart
00179 !
00180 #ifdef PRISM_ASSERTION
00181       if (n /= n_send) then
00182          print *, 'n, n_send', n, n_send
00183          call psmile_assert ( __FILE__, __LINE__, &
00184                  'All indices were NOT generated for srcloc_ind')
00185       endif
00186 #endif
00187 !
00188 !===> All done
00189 !
00190 #ifdef VERBOSE
00191       print 9980, trim(ch_id), ierror
00192 
00193       call psmile_flushstd
00194 #endif /* VERBOSE */
00195 !
00196 !  Formats:
00197 !
00198 #ifdef VERBOSE
00199 
00200 9990 format (1x, a, ': psmile_get_face_ind_3d:')
00201 9980 format (1x, a, ': psmile_get_face_ind_3d: eof ierror =', i3)
00202 
00203 #endif /* VERBOSE */
00204 
00205 #ifdef DEBUG
00206 #endif
00207 
00208       end subroutine PSMILe_Get_face_ind_3d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1