psmile_store_source_locs_3d.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_3d (found, loc, ibeg, len, &
00007                       send_info, nloc, opt,                          &
00008                       nprev, nadd, ierror)
00009 !
00010 ! !USES:
00011 !
00012       use PRISM_constants
00013 !
00014       use PSMILe, dummy_interface => PSMILe_store_source_locs_3d
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,k".
00031 !     Level number = nlev+1: Never found (input value)
00032 !
00033       Integer, Intent (In)            :: loc (ndim_3d, len)
00034 
00035 !     Locations found for coordinates to be searched
00036 
00037       Integer, Intent (In)            :: nloc
00038 !
00039 !     Total 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)            :: nprev
00051 !
00052 !     Number of locations already stored in
00053 !     send_info%srclocs(ialloc, ipart)%vector
00054 !
00055 ! !INPUT/OUTPUT PARAMETERS:
00056 !
00057       Type(Send_information), Intent(InOut) :: send_info
00058 !
00059 !     Send information to be generated
00060 !
00061 ! !OUTPUT PARAMETERS:
00062 !
00063       Integer, Intent (Out)           :: nadd
00064 !
00065 !     Number of locations added in send_info%srclocs(ialloc, ipart)%vector
00066 !
00067       integer, Intent (Out)           :: ierror
00068 
00069 !     Returns the error code of PSMILe_Store_source_locs_3d;
00070 !             ierror = 0 : No error
00071 !             ierror > 0 : Severe error
00072 !
00073 ! !LOCAL VARIABLES
00074 !
00075       Integer, Parameter              :: ialloc = 1
00076       Integer, Parameter              :: ipart  = 1
00077 !
00078       Integer                         :: i, n
00079 !
00080       Integer, Pointer                :: srcloc (:)
00081 !
00082 !     ... for error handling
00083 !
00084       Integer, Parameter              :: nerrp = 2
00085       Integer                         :: ierrp (nerrp)
00086 !
00087 ! !DESCRIPTION:
00088 !
00089 ! Subroutine "PSMILe_store_source_locs_3d" stores the data on the
00090 ! source locations found for the method (grid) and the
00091 ! subgrid coords sent by sending (destination) process.
00092 !
00093 !
00094 ! !REVISION HISTORY:
00095 !
00096 !   Date      Programmer   Description
00097 ! ----------  ----------   -----------
00098 ! 03.07.03    H. Ritzdorf  created
00099 !
00100 !EOP
00101 !----------------------------------------------------------------------
00102 !
00103 !  $Id: psmile_store_source_locs_3d.F90 2325 2010-04-21 15:00:07Z valcke $
00104 !  $Author: valcke $
00105 !
00106    Character(len=len_cvs_string), save :: mycvs = 
00107        '$Id: psmile_store_source_locs_3d.F90 2325 2010-04-21 15:00:07Z valcke $'
00108 !
00109 !----------------------------------------------------------------------
00110 !
00111 !  Initialization
00112 !
00113 #ifdef VERBOSE
00114       print 9990, trim(ch_id), opt, nloc, ibeg, len, nprev
00115 
00116       call psmile_flushstd
00117 #endif /* VERBOSE */
00118 !
00119       ierror = 0
00120 !
00121       if (nprev == 0 .and. &
00122           .not. Associated (send_info%srclocs(ialloc, ipart)%vector) ) then
00123 
00124 !        Allocate (send_info%srcloc(1:ndim_3d, 1:nloc), STAT = ierror)
00125          Allocate (send_info%srclocs(ialloc, ipart)%vector(ndim_3d*nloc), &
00126                    STAT = ierror)
00127          if ( ierror > 0 ) then
00128             ierrp (1) = ierror
00129             ierrp (2) = nloc * ndim_3d
00130 
00131             ierror = PRISM_Error_Alloc
00132             call psmile_error ( ierror, 'send_info%srclocs(1)%vector', &
00133                                 ierrp, 2, __FILE__, __LINE__ )
00134             return
00135 
00136          endif
00137 !
00138 !        Allocate dummy array (for 1 dummy area) in order to avoid
00139 !        problems with some run-time systems.
00140 !
00141          Allocate (send_info%srcars(ialloc, ipart)%vector(ndim_3d*2), &
00142                    STAT = ierror)
00143          if ( ierror > 0 ) then
00144             ierrp (1) = ierror
00145             ierrp (2) = ndim_3d * 2
00146 
00147             ierror = PRISM_Error_Alloc
00148             call psmile_error ( ierror, 'send_info%srcars(1)%vector', &
00149                                 ierrp, 2, __FILE__, __LINE__ )
00150             return
00151 
00152          endif
00153 
00154 !
00155 ! ... Store locations in a simple list
00156 !     ??? Schoener waere eine clusterung in Regions
00157 !     ??? Dann muessen aber die Destionations mit den geclusterten
00158 !         Daten erzeugt werden, um Problem in der Reihenfolge zu verhindern.
00159 !
00160          send_info%npoints(ialloc, ipart) = 0
00161          send_info%nars   (ialloc, ipart) = 0
00162       end if
00163 !
00164       srcloc => send_info%srclocs(ialloc, ipart)%vector
00165 !
00166       n = nprev * ndim_3d
00167 
00168       if (opt == 0) then
00169 !cdir vector
00170           do i = ibeg, len
00171           if (abs(found (i)) == 1) then
00172 
00173              srcloc (n+1) = loc (1, i)
00174              srcloc (n+2) = loc (2, i)
00175              srcloc (n+3) = loc (3, i)
00176 
00177              n = n + ndim_3d
00178 
00179 #ifdef DEBUGX
00180              if (n/ndim_3d == 1) then
00181                 print *, 'store_source_locs_3d: n, i, srcloc', &
00182                          n/ndim_3d, i, srcloc (n-2:n), loc (:, i), found (i)
00183              endif
00184 #endif
00185           endif
00186           end do
00187 !
00188       else
00189 !
00190 !cdir vector
00191           do i = ibeg, len
00192           if (found (i) == opt) then
00193              srcloc (n+1) = loc (1, i)
00194              srcloc (n+2) = loc (2, i)
00195              srcloc (n+3) = loc (3, i)
00196 
00197              n = n + ndim_3d
00198 
00199 #ifdef DEBUGX
00200              if (n/ndim_3d == 3441) then
00201                 print *, 'store_source_locs_3d: n',  n/ndim_3d, ', i', i, &
00202                          'srcloc', srcloc (n-2:n), &
00203                          'loc', loc (:, i), ', found', found (i)
00204              endif
00205 #endif
00206           endif
00207           end do
00208       endif
00209 
00210 #ifdef PRISM_ASSERTION
00211       if (nloc*ndim_3d < n) then
00212          print *, "nloc, n", nloc, n/ndim_3d
00213          call psmile_assert ( __FILE__, __LINE__, "nloc < n")
00214       endif
00215 #endif
00216 !
00217 !===> All done
00218 !
00219       nadd = n/ndim_3d - nprev
00220 !
00221       send_info%npoints(ialloc, ipart) = send_info%npoints(ialloc, ipart) + nadd
00222 !     send_info%nars   (ialloc, ipart) = send_info%nars   (ialloc, ipart) +
00223 !
00224 #ifdef VERBOSE
00225       print 9980, trim(ch_id), ierror
00226 
00227       call psmile_flushstd
00228 #endif /* VERBOSE */
00229 !
00230 !  Formats:
00231 !
00232 #ifdef VERBOSE
00233 
00234 9990 format (1x, a, ': psmile_store_source_locs_3d: opt ', i2, &
00235              '; nloc, ibeg, len, nprev ', 4i8)
00236 9980 format (1x, a, ': psmile_store_source_locs_3d: eof ierror =', i3)
00237 
00238 #endif /* VERBOSE */
00239 
00240       end subroutine PSMILe_store_source_locs_3d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1