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

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1