psmile_store_dest_locs_21d.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2007-2010, NEC Europe Ltd., London, UK.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006       subroutine psmile_store_dest_locs_21d (found, range,  control,  &
00007                                              foundz, &
00008                                              send_info, nloc, opt,           &
00009                                              nprev, nadd, ierror)
00010 !
00011 ! !USES:
00012 !
00013       use PRISM_constants
00014 !
00015       use PSMILe, dummy_interface => PSMILe_Store_dest_locs_21d
00016 
00017       implicit none
00018 !
00019 ! !INPUT PARAMETERS:
00020 !
00021       Integer, Intent (In)            :: range (2, ndim_3d)
00022 
00023 !     Dimension of loc and found
00024 
00025       Integer, Intent (In)            :: control  (2, ndim_3d)
00026 
00027 !     Index range found
00028 
00029       Integer, Intent (In)            :: found  (range(1,1):range(2,1), 
00030                                                  range(1,2):range(2,2))
00031       Integer, Intent (In)            :: foundz (range(1,3):range(2,3))
00032 
00033 !     Finest level number on which a grid cell was found for point i,j,k.
00034 !     Level number = nlev+1: Never found (input value)
00035 !
00036       Integer, Intent (In)            :: nloc
00037 !
00038 !     Number of locations to be stored
00039 
00040       Integer, Intent (In)            :: opt
00041 !
00042 !     Option:
00043 !     = -1 : Store values which have to be sent to the coupler;
00044 !            i.e. store locations with found (i,j,k) = -1
00045 !     =  1 : Store values which have to be sent to the application process;
00046 !            i.e. store locations with found (i,j,k) =  1
00047 !     =  0 : Store values with abs(found (i,j,k)) = 1.
00048 !
00049       Integer, Intent (In)            :: nprev
00050 !
00051 !     Number of locations already stored in send_info%dstijk
00052 !
00053 !
00054 ! !OUTPUT PARAMETERS:
00055 !
00056       Type(Send_information), Intent(Inout) :: send_info
00057 !
00058 !     Send information to be generated
00059 !
00060       Integer, Intent (Out)           :: nadd
00061 !
00062 !     Number of locations added in send_info%dstijk
00063 !
00064       integer, Intent (Out)           :: ierror
00065 
00066 !     Returns the error code of PSMILe_Store_dest_locs_21d;
00067 !             ierror = 0 : No error
00068 !             ierror > 0 : Severe error
00069 !
00070 ! !LOCAL VARIABLES
00071 !
00072       Integer, Parameter              :: val_direct  =  1
00073       Integer, Parameter              :: val_coupler = -1
00074       Integer, Parameter              :: val_both    =  0
00075       Integer, Parameter              :: val_abs     =  1
00076 !
00077 !     ... Loop variables
00078 !
00079       Integer                         :: i, j, k
00080 !
00081 !     ... Counters for number of points
00082 !
00083       Integer                         :: n, noldk, nrep
00084       Logical                         :: anyk
00085 !
00086 !     ... for error handling
00087 !
00088       Integer, Parameter              :: nerrp = 2
00089       Integer                         :: ierrp (nerrp)
00090 !
00091 ! !DESCRIPTION:
00092 !
00093 ! Subroutine "PSMILe_Store_dest_locs_21d" stores the data on the
00094 ! the subgrid coords sent by sending (destination) process which were found
00095 ! and which will be sent by the actual process.
00096 !
00097 ! Subroutine "PSMILe_Store_dest_locs_21d" is designed for grids of
00098 ! of type "PRISM_Irrlonlat_regvrt".
00099 !
00100 ! !REVISION HISTORY:
00101 !
00102 !   Date      Programmer   Description
00103 ! ----------  ----------   -----------
00104 ! 03.07.21    H. Ritzdorf  created
00105 !
00106 !EOP
00107 !----------------------------------------------------------------------
00108 !
00109 !  $Id: psmile_store_dest_locs_21d.F90 2687 2010-10-28 15:15:52Z coquart $
00110 !  $Author: coquart $
00111 !
00112    Character(len=len_cvs_string), save :: mycvs = 
00113        '$Id: psmile_store_dest_locs_21d.F90 2687 2010-10-28 15:15:52Z coquart $'
00114 !
00115 !----------------------------------------------------------------------
00116 !
00117 !  Initialization
00118 !
00119 #ifdef VERBOSE
00120       print 9990, trim(ch_id), opt, control
00121 
00122       call psmile_flushstd
00123 #endif /* VERBOSE */
00124 !
00125       ierror = 0
00126 !
00127       if (nprev == 0 .and. .not. Associated (send_info%dstijk) ) then
00128 
00129          Allocate (send_info%dstijk(1:ndim_3d, 1:nloc), STAT = ierror)
00130          if ( ierror > 0 ) then
00131             ierrp (1) = ierror
00132             ierrp (2) = nloc * ndim_3d
00133 
00134             ierror = PRISM_Error_Alloc
00135 
00136             call psmile_error ( ierror, 'send_info%dstijk', &
00137                                 ierrp, 2, __FILE__, __LINE__ )
00138             return
00139          endif
00140       endif
00141 !
00142 ! ... Store locations in a simple list
00143 ! ... Schoener waere eine clusterung in Regions
00144 
00145       n = nprev
00146       nrep = 0
00147 
00148       if (opt == val_both) then
00149 !
00150 !     ... Find points which have to be sent
00151 !     ... Generate data for 2d hyperplane (first 2 indices)
00152 !
00153          do j = control (1, 2), control (2, 2)
00154 !cdir vector
00155             do i = control (1, 1), control (2, 1)
00156             if (abs(found (i,j)) == val_abs) then
00157                n = n + 1
00158 
00159                send_info%dstijk (1, n) = i
00160                send_info%dstijk (2, n) = j
00161             endif
00162             end do
00163          end do
00164 !
00165          nadd = n - nprev
00166 !
00167 !     ... Generate data for 3rd (regular) index
00168 !
00169 !cdir vector
00170           do k = control (1, 3), control (2, 3)
00171             if (abs(foundz (k)) == val_abs) then
00172                send_info%dstijk (3, nprev+nrep*nadd+1:nprev+(nrep+1)*nadd) = k
00173                nrep = nrep + 1
00174             endif
00175           end do
00176 !
00177 !     ...  Replicate data for 2d Hyperplane; i.e.
00178 !          copy first 2 indices for each 3rd index
00179 !
00180          do i = 1, nrep-1
00181             send_info%dstijk (1:2, nprev+i*nadd+1:nprev+(i+1)*nadd) = &
00182             send_info%dstijk (1:2, nprev       +1:nprev+nadd)
00183          end do
00184 !
00185          nadd = nadd * nrep
00186 !
00187       else if (opt == val_direct) then
00188 !
00189 !     ... Find points which have to be sent directly to the destination process
00190 !     ... Generate data for 2d hyperplane (first 2 indices)
00191 !
00192          do j = control (1, 2), control (2, 2)
00193 !cdir vector
00194             do i = control (1, 1), control (2, 1)
00195             if (found (i,j) == val_direct) then
00196                n = n + 1
00197 
00198                send_info%dstijk (1, n) = i
00199                send_info%dstijk (2, n) = j
00200             endif
00201             end do
00202          end do
00203 !
00204          nadd = n - nprev
00205 !
00206 !     ... Generate data for 3rd (regular) index
00207 !
00208 !cdir vector
00209           do k = control (1, 3), control (2, 3)
00210             if (foundz (k) == val_direct) then
00211                send_info%dstijk (3, nprev+nrep*nadd+1:nprev+(nrep+1)*nadd) = k
00212                nrep = nrep + 1
00213             endif
00214           end do
00215 !
00216 !     ...  Replicate data for 2d Hyperplane; i.e.
00217 !          copy first 2 indices for each 3rd index
00218 !
00219          do i = 1, nrep-1
00220             send_info%dstijk (1:2, nprev+i*nadd+1:nprev+(i+1)*nadd) = &
00221             send_info%dstijk (1:2, nprev       +1:nprev+nadd)
00222          end do
00223 !
00224          nadd = nadd * nrep
00225 !
00226       else if (opt == val_coupler) then
00227 !
00228 !     ... Find points which have to be sent to the coupler
00229 !
00230 !cdir vector
00231          do k = control (1, 3), control (2, 3)
00232             if (abs(foundz(k)) /= val_abs) continue
00233 
00234             anyk = foundz(k) == val_coupler
00235             noldk = n
00236 !
00237 !     ... Generate data for 2d hyperplane (first 2 indices)
00238 !
00239             do j = control (1, 2), control (2, 2)
00240 !cdir vector
00241                do i = control (1, 1), control (2, 1)
00242                   if (found(i,j) == val_coupler .or. &
00243                      (anyk .and. found(i,j) == val_direct)) then
00244                      n = n + 1
00245 
00246                      send_info%dstijk (1, n) = i
00247                      send_info%dstijk (2, n) = j
00248                   endif
00249                end do
00250             end do
00251 !
00252 !     ... Generate data for 3rd (regular) index
00253 !
00254             send_info%dstijk (3, noldk+1:n) = k
00255          end do
00256 !
00257          nadd = n - nprev
00258 !
00259       endif
00260 !
00261 #ifdef PRISM_ASSERTION
00262 !
00263 !     Control correct dimension
00264 !
00265       if (nloc < nprev+nadd) then
00266          print *, "nloc, nctl", nloc, nprev+nadd, nrep
00267          call psmile_assert ( __FILE__, __LINE__, "nloc < nprev+nadd")
00268       endif
00269 #endif
00270 !
00271 !===> All done
00272 !
00273 #ifdef VERBOSE
00274       print 9980, trim(ch_id), ierror
00275 
00276       call psmile_flushstd
00277 #endif /* VERBOSE */
00278 !
00279 !  Formats:
00280 !
00281 #ifdef VERBOSE
00282 
00283 9990 format (1x, a, ': psmile_store_dest_locs_21d: opt ', i2, &
00284              '; control ', 6i6)
00285 9980 format (1x, a, ': psmile_store_dest_locs_21d: eof ierror =', i3)
00286 
00287 #endif /* VERBOSE */
00288 
00289       end subroutine PSMILe_Store_dest_locs_21d

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1