psmile_store_source_locs_1d.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       subroutine psmile_store_source_locs_1d (found,  loc,  ibeg,  len,  &
00007                                               send_info, nloc, opt,      &
00008                                               ialloc, ipart, nprev, nadd, ierror)
00009 !
00010 ! !USES:
00011 !
00012       use PRISM_constants
00013 !
00014       use PSMILe, dummy_interface => PSMILe_store_source_locs_1d
00015 
00016       implicit none
00017 !
00018 ! !INPUT PARAMETERS:
00019 !
00020       Integer, Intent (In)            :: ibeg
00021 !
00022 !     Start of search
00023 !
00024       Integer, Intent (In)            :: len
00025 !
00026 !     Dimension of vectors found and loc
00027 !
00028       Integer, Intent (In)            :: found  (len)
00029 
00030 !     Finest level number on which a grid cell was found for point i,j.
00031 !     Level number = nlev+1: Never found (input value)
00032 !
00033       Integer, Intent (In)            :: loc  (len)
00034 
00035 !     Locations found for coordinates to be searched
00036 
00037       Integer, Intent (In)            :: nloc
00038 !
00039 !     Number of locations to be stored
00040 
00041       Integer, Intent (In)            :: opt
00042 !
00043 !     Option:
00044 !     =  1 : Store values which have to be sent to the coupler;
00045 !            i.e. store locations with found (i,j,k) =  1
00046 !     = -1 : Store values which have to be sent to the application process;
00047 !            i.e. store locations with found (i,j,k) = -1
00048 !     =  0 : Store values with abs(found (i,j,k)) = 1.
00049 
00050       Integer, Intent (In)            :: ialloc
00051 !
00052 !     Index of source locations to be created in "send_info".
00053 
00054       Integer, Intent (In)            :: ipart
00055 !
00056 !     Number of the partition
00057 
00058       Integer, Intent (In)            :: nprev
00059 !
00060 !     Number of locations already stored in send_info%srclocs(ialloc, ipart)%vector
00061 !
00062 ! !OUTPUT PARAMETERS:
00063 !
00064       Type(Send_information), Intent(InOut) :: send_info
00065 !
00066 !     Send information to be generated
00067 !
00068       Integer, Intent (Out)           :: nadd
00069 !
00070 !     Number of locations added in send_info%srclocs(ialloc, ipart)%vector
00071 !
00072       integer, Intent (Out)           :: ierror
00073 
00074 !     Returns the error code of PSMILe_Store_source_locs_1d;
00075 !             ierror = 0 : No error
00076 !             ierror > 0 : Severe error
00077 !
00078 ! !LOCAL VARIABLES
00079 !
00080       Integer                         :: i, n
00081 !
00082       Integer, Pointer                :: srcloc (:)
00083 !
00084 !     ... for error handling
00085 !
00086       Integer, Parameter              :: nerrp = 2
00087       Integer                         :: ierrp (nerrp)
00088 !
00089 ! !DESCRIPTION:
00090 !
00091 ! Subroutine "PSMILe_store_source_locs_1d" stores the data on the
00092 ! source locations found for the method (grid) and the
00093 ! subgrid coords sent by sending (destination) process.
00094 !
00095 ! Subroutine "PSMILe_Store_dest_locs_1d" is designed for (1d regular sub-)
00096 ! grids of of type "PRISM_Irrlonlat_regvrt".
00097 !
00098 ! !REVISION HISTORY:
00099 !
00100 !   Date      Programmer   Description
00101 ! ----------  ----------   -----------
00102 ! 03.07.03    H. Ritzdorf  created
00103 !
00104 !EOP
00105 !----------------------------------------------------------------------
00106 !
00107 !  $Id: psmile_store_source_locs_1d.F90 2325 2010-04-21 15:00:07Z valcke $
00108 !  $Author: valcke $
00109 !
00110    Character(len=len_cvs_string), save :: mycvs = 
00111        '$Id: psmile_store_source_locs_1d.F90 2325 2010-04-21 15:00:07Z valcke $'
00112 !
00113 !----------------------------------------------------------------------
00114 !
00115 !  Initialization
00116 !
00117 #ifdef VERBOSE
00118       print 9990, trim(ch_id), opt, ibeg, len
00119 
00120       call psmile_flushstd
00121 #endif /* VERBOSE */
00122 !
00123       ierror = 0
00124 !
00125       if (nprev == 0 .and. &
00126           .not. Associated (send_info%srclocs(ialloc, ipart)%vector) ) then
00127 
00128          Allocate (send_info%srclocs(ialloc, ipart)%vector(nloc), STAT = ierror)
00129          if ( ierror > 0 ) then
00130             ierrp (1) = ierror
00131             ierrp (2) = nloc
00132             call psmile_error ( PRISM_Error_Alloc, 'send_info%srclocs()%vector', &
00133                                 ierrp, 2, __FILE__, __LINE__ )
00134             return
00135          endif
00136 
00137       endif
00138 !
00139 ! ... Store locations in a simple list
00140 !     ??? Schoener waere eine clusterung in Regions
00141 !     ??? Dann muessen aber die Destionations mit den geclusterten
00142 !         Daten erzeugt werden, um Problem in der Reihenfolge zu verhindern.
00143 !
00144       srcloc => send_info%srclocs(ialloc, ipart)%vector
00145 !
00146       n = nprev
00147 
00148       if (opt == 0) then
00149 !cdir vector
00150           do i = ibeg, len
00151           if (abs(found (i)) == 1) then
00152              srcloc (n+1) = loc (i)
00153              n = n + 1
00154           endif
00155           end do
00156 !
00157       else
00158 !
00159 !cdir vector
00160           do i = ibeg, len
00161           if (found (i) == opt) then
00162              srcloc (n+1) = loc (i)
00163              n = n + 1
00164           endif
00165           end do
00166       endif
00167 
00168 #ifdef PRISM_ASSERTION
00169       if (nloc < n) then
00170          print *, "nloc, n", nloc, n
00171          call psmile_assert ( __FILE__, __LINE__, "nloc < n")
00172       endif
00173 #endif
00174 !
00175 !===> All done
00176 !
00177       nadd = n - nprev
00178 !
00179       send_info%npoints(ialloc, ipart) = send_info%npoints(ialloc, ipart) + nadd
00180 !     send_info%nars   (ialloc, ipart) = send_info%nars   (ialloc, ipart) +
00181 
00182 #ifdef DEBUGX
00183       n = nprev
00184       if (opt == 0) then
00185          do i = ibeg, len
00186          if (abs(found (i)) == 1) then
00187             n = n + 1
00188             print *, i, loc(i), n, srcloc(n)
00189          endif
00190          end do
00191       else
00192          do i = ibeg, len
00193          if (found (i) == opt) then
00194             n = n + 1
00195             print *, i, loc(i), n, srcloc(n)
00196          endif
00197          end do
00198       endif
00199 #endif
00200 !
00201 #ifdef VERBOSE
00202       print 9980, trim(ch_id), ierror
00203 
00204       call psmile_flushstd
00205 #endif /* VERBOSE */
00206 !
00207 !  Formats:
00208 !
00209 #ifdef VERBOSE
00210 
00211 9990 format (1x, a, ': psmile_store_source_locs_1d: opt ', i2, &
00212              '; ibeg, len ', 2i6)
00213 9980 format (1x, a, ': psmile_store_source_locs_1d: eof ierror =', i3)
00214 
00215 #endif /* VERBOSE */
00216 
00217       end subroutine PSMILe_store_source_locs_1d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1