psmile_get_face_ind_reg.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_face_ind_reg
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_get_face_ind_reg (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_reg
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_reg;
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, nprev
00079       Integer                         :: leni, lenij
00080 !
00081       Integer, Pointer                :: indices_req (:)
00082       Integer, Pointer                :: len_req (:)
00083       Integer, Pointer                :: srcloc_i (:), srcloc_j (:)
00084       Integer, Pointer                :: srcloc_k (:)
00085 !
00086 !     ... for indices
00087 !
00088       Integer                         :: ind (ndim_3d)
00089 !
00090 !     ... for error handling
00091 !
00092 !     Integer, Parameter              :: nerrp = 3
00093 !     Integer                         :: ierrp (nerrp)
00094 !
00095 ! !DESCRIPTION:
00096 !
00097 ! Subroutine "PSMILe_Get_face_ind_reg" computes the 3d indices of
00098 ! the locations to be sent to a "coinciding" (neighboured)
00099 ! process of the same component.
00100 !
00101 !
00102 ! !REVISION HISTORY:
00103 !
00104 !   Date      Programmer   Description
00105 ! ----------  ----------   -----------
00106 ! 02.02.05    H. Ritzdorf  created
00107 !
00108 !EOP
00109 !----------------------------------------------------------------------
00110 !
00111 !  $Id: psmile_get_face_ind_reg.F90 2082 2009-10-21 13:31:19Z hanke $
00112 !  $Author: hanke $
00113 !
00114    Character(len=len_cvs_string), save :: mycvs = 
00115        '$Id: psmile_get_face_ind_reg.F90 2082 2009-10-21 13:31:19Z hanke $'
00116 !
00117 !----------------------------------------------------------------------
00118 !
00119 !  Initialization
00120 !
00121 #ifdef VERBOSE
00122       print 9990, trim(ch_id)
00123 
00124       call psmile_flushstd
00125 #endif /* VERBOSE */
00126 !
00127       ierror = 0
00128 !
00129       len_req => extra_search%len_req
00130 !
00131 #ifdef PRISM_ASSERTION
00132       if (SUM(len_req(:)) > nreq) then
00133          print *, 'nreq, sum(len_req)', nreq, SUM(len_req(:))
00134          call psmile_assert ( __FILE__, __LINE__, &
00135                  'nreq is too small (inconsistent data)')
00136       endif
00137 !
00138       if (send_info%nvec /= ndim_3d) then
00139          print *, 'nvec', send_info%nvec
00140          call psmile_assert ( __FILE__, __LINE__, &
00141                  'Routine is designed for regular case (3 1d vectors)')
00142       endif
00143 #endif
00144 !
00145 !  For all parts "ipart"
00146 !
00147 !  indices_req = Indices of the extra points in entire list of all points
00148 !                to be searched (i.e. for all parts "ipart")
00149 !
00150       n = 0
00151       nprev = 0
00152       nextra_prev = 0
00153 ! 
00154          do ipart = 1, search%npart
00155 !
00156          if (len_req (ipart) > 0) then
00157             indices_req => extra_search%indices_req(ipart)%vector
00158             srcloc_i => send_info%srclocs(1,ipart)%vector
00159             srcloc_j => send_info%srclocs(2,ipart)%vector
00160             srcloc_k => send_info%srclocs(3,ipart)%vector
00161 !
00162             leni  = send_info%npoints(1,ipart)
00163             lenij = send_info%npoints(2,ipart) * leni
00164 !
00165                do j = 1, len_req (ipart)
00166                if (send_mask(nextra_prev+j)) then
00167 !
00168 !  ... compute 3d index out of 1-dimensional index "i"
00169 !      (srclocs are stored as 3 1-dimensional vectors
00170 !       regulary in lon, lat and z-direction)
00171 !
00172                   n = n + 1
00173 !
00174                   i = indices_req(j) - nprev - 1
00175 !
00176                   ind (3) =  i / lenij
00177                   ind (2) = (i-ind(3)*lenij) / leni
00178                   ind (1) = mod (i, leni)
00179 !
00180 #ifdef PRISM_ASSERTION
00181                   if (ind(1) < 0 .or. ind(1) >= send_info%npoints(1,ipart) .or. &
00182                       ind(2) < 0 .or. ind(2) >= send_info%npoints(2,ipart) .or. &
00183                       ind(3) < 0 .or. ind(3) >= send_info%npoints(3,ipart) ) then
00184                      print *, "j, i, ind", j, i, ind
00185                      call psmile_assert (__FILE__, __LINE__, &
00186                         "Incorrect index generated");
00187                   endif
00188 #endif
00189 !
00190                   srcloc_ind (1, n) = srcloc_i(ind(1)+1)
00191                   srcloc_ind (2, n) = srcloc_j(ind(2)+1)
00192                   srcloc_ind (3, n) = srcloc_k(ind(3)+1)
00193 print *, 'get_face_ind: n', n, srcloc_ind (:, n)
00194                endif
00195                end do ! j
00196 !
00197             nextra_prev = nextra_prev + len_req (ipart)
00198          endif
00199 !
00200          nprev = nprev + len_cpl (ipart)
00201          end do ! ipart
00202 !
00203 #ifdef PRISM_ASSERTION
00204       if (n /= n_send) then
00205          print *, 'n, n_send', n, n_send
00206          call psmile_assert ( __FILE__, __LINE__, &
00207                  'All indices were NOT generated for srcloc_ind')
00208       endif
00209 #endif
00210 !
00211 !===> All done
00212 !
00213 #ifdef VERBOSE
00214       print 9980, trim(ch_id), ierror
00215 
00216       call psmile_flushstd
00217 #endif /* VERBOSE */
00218 !
00219 !  Formats:
00220 !
00221 #ifdef VERBOSE
00222 
00223 9990 format (1x, a, ': psmile_get_face_ind_reg:')
00224 9980 format (1x, a, ': psmile_get_face_ind_reg: eof ierror =', i3)
00225 
00226 #endif /* VERBOSE */
00227 
00228 #ifdef DEBUG
00229 #endif
00230 
00231       end subroutine PSMILe_Get_face_ind_reg

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1