psmile_mg_method_1d_dble.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_MG_method_1d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_method_1d_dble (                &
00012                            comp_info, nlev,                &
00013                            found, loc, range,              &
00014                            coords1,                        &
00015                            shape, control,                 &
00016                            method_id,                      &
00017                            x_coords,                       &
00018                            coords_shape,                   &
00019                            grid_valid_shape, cyclic,       &
00020                            chmin, chmax,                   &
00021                            tol, ierror)
00022 !
00023 ! !USES:
00024 !
00025       use PRISM_constants
00026 !
00027       use PSMILe, dummy_interface => PSMILe_mg_method_1d_dble
00028 
00029       implicit none
00030 !
00031 ! !INPUT PARAMETERS:
00032 !
00033       Type (Enddef_comp), Intent (In) :: comp_info
00034 
00035 !     Info on the component in which the donor
00036 !     should be searched.
00037 
00038       Integer, Intent (In)            :: nlev
00039 
00040 !     Total number of levels
00041 
00042       Integer, Intent (In)            :: range (2, ndim_3d)
00043 
00044 !     Dimension of loc and found
00045 
00046       Integer, Intent (In)            :: shape (2, ndim_3d)
00047 
00048 !     Dimension of coordinate arrays
00049 
00050       Integer, Intent (InOut)         :: found (range(1,1):range(2,1), 
00051                                                 range(1,2):range(2,2), 
00052                                                 range(1,3):range(2,3))
00053 
00054 !     Finest level number on which a grid cell was found for point i,j,k.
00055 !     Input :
00056 !     Level number < -1: Point was not found and
00057 !                        and last level number was (-found(i,j,k))
00058 !     Level number = -(nlev+1): Never found
00059 !     Output :
00060 !     found(i,j,k) = +1: Point (i,j,k) is located on a point
00061 !                        of the method grid.
00062 !     found(i,j,k) = -1: Point (i,j,k) is located in a cell
00063 !                        of the method grid.
00064 !     Otherwise, not contained in the method grid.
00065 
00066       Integer, Intent (InOut)         :: loc   (range(1,1):range(2,1), 
00067                                                 range(1,2):range(2,2), 
00068                                                 range(1,3):range(2,3))
00069 
00070 !     Indices of the grid cell in which the point was found.
00071 !     The indices are relative to "shape".
00072 
00073       Double Precision, Intent (In)   :: coords1 (shape(1,1):shape(2,1), 
00074                                                   shape(1,2):shape(2,2), 
00075                                                   shape(1,3):shape(2,3))
00076 
00077 !     Coordinates to be searched
00078 
00079       Integer, Intent (In)            :: control (2, ndim_3d)
00080 
00081 !     Input  value: index range in "coords" to be searched
00082 !     Output value: index range in "coords" for which a cell was found
00083 
00084       Integer, Intent (In)            :: method_id
00085 
00086 !     Method id
00087 
00088       Integer, Intent (In)            :: coords_shape (2)
00089 
00090 !     Dimension of coordinates (method) x_coords, ...
00091 
00092       Double Precision, Intent (In)   :: x_coords(coords_shape(1):coords_shape(2))
00093 
00094 !     Coordinates of the method
00095 
00096       Integer,          Intent (In)   :: grid_valid_shape (2)
00097 
00098 !     Specifies the valid block shape for the "1d"-dimensional block
00099 !     Dimension of chmin, chmax and midp
00100 
00101       Logical,          Intent (In)   :: cyclic
00102 
00103 !     Is this a (locally) cyclic coordinate
00104 
00105       Double Precision, Intent (In)   :: chmin (grid_valid_shape(1)-1:
00106                                                 grid_valid_shape(2))
00107 !
00108 !     Minimum value of bounding box of coordinate.
00109 ! chmin (grid_valid_shape(1)+i) contains the minimal value of bounding box of
00110 !                   cell (x_coords(grid_valid_shape(1)+i)  :
00111 !                         x_coords(grid_valid_shape(1)+i+1))
00112 !
00113 ! The first box
00114 ! chmin (grid_valid_shape(1)-1) contains the minimal value of
00115 !                                        the bounding box from
00116 !                    the virtual cell (x_coords(grid_valid_shape(1)):
00117 !                                               bounding box of
00118 !                                               corner control volume)
00119 ! The last box
00120 ! chmin (grid_valid_shape(2)) contains the minimal value of
00121 !                                      the bounding box from
00122 !                    the virtual cell (x_coords(grid_valid_shape(2)):
00123 !                                               bounding box of
00124 !                                               corner control volume)
00125 !
00126       Double Precision, Intent (In)   :: chmax (grid_valid_shape(1)-1:
00127                                                 grid_valid_shape(2))
00128 !
00129 !     Maximum value of bounding box of coordinate.
00130 !
00131       Double Precision, Intent (In)   :: tol
00132 !
00133 !     Absolute tolerance for search of "identical" points
00134 !     TOL >= 0.0
00135 !
00136 ! !OUTPUT PARAMETERS:
00137 !
00138       integer, Intent (Out)           :: ierror
00139 
00140 !     Returns the error code of PSMILE_mg_method_1d_dble;
00141 !             ierror = 0 : No error
00142 !             ierror > 0 : Severe error
00143 !
00144 ! !DEFINED PARAMETERS:
00145 !
00146 !  NREFD = Number of method coordinates to be controlled.
00147 !  val_direct  = Code for locations which should to be directly transferred
00148 !                to the destination process.
00149 !  val_coupler = Code for locations which should to be transferred
00150 !                to the coupler process.
00151 
00152       Integer, Parameter              :: nrefd = 3
00153 !
00154       Integer, Parameter              :: val_direct  =  1
00155       Integer, Parameter              :: val_coupler = -1
00156 !
00157 !     ... for levels
00158 !
00159       Integer, Parameter              :: lev = 1
00160 !
00161 ! !LOCAL VARIABLES
00162 !
00163 !     ... for locations searched
00164 !
00165       Integer                         :: i, j, k
00166       Integer                         :: ibegl, iendl, jpart, n
00167       Integer                         :: ic
00168 !
00169 !     ... For locations controlled
00170 !
00171       Double Precision                :: dist (nrefd)
00172       Integer                         :: irange (2)
00173 !
00174       Integer                         :: nref
00175       Integer                         :: nmin (1)
00176 !
00177 !     ... for search
00178 !
00179       Integer                         :: ii
00180 !
00181 #ifdef DEBUG
00182 !     ... for locations searched
00183 !
00184       Integer                         :: nfnd (0:3)
00185 #endif /* DEBUG */
00186 !
00187 ! !DESCRIPTION:
00188 !
00189 ! Subroutine "PSMILe_mg_method_1d_dble" searches the corresponding points
00190 ! of the method (grid) for the subgrid coords by the sending process.
00191 !
00192 ! TODO: Besser machen !!
00193 !
00194 ! !REVISION HISTORY:
00195 !
00196 !   Date      Programmer   Description
00197 ! ----------  ----------   -----------
00198 ! 03.07.21    H. Ritzdorf  created
00199 !
00200 !EOP
00201 !----------------------------------------------------------------------
00202 !
00203 !  $Id: psmile_mg_method_1d_dble.F90 2325 2010-04-21 15:00:07Z valcke $
00204 !  $Author: valcke $
00205 !
00206    Character(len=len_cvs_string), save :: mycvs = 
00207        '$Id: psmile_mg_method_1d_dble.F90 2325 2010-04-21 15:00:07Z valcke $'
00208 !----------------------------------------------------------------------
00209 !
00210 !  Initialization
00211 !
00212       ierror = 0
00213 
00214 #ifdef VERBOSE
00215       print 9990, trim(ch_id), lev, control
00216 
00217       call psmile_flushstd
00218 #endif /* VERBOSE */
00219 !
00220 #ifdef PRISM_ASSERTION
00221       if (tol < 0.0) then
00222          call psmile_assert (__FILE__, __LINE__, &
00223                              "tol must be >= 0.0")
00224       endif
00225 
00226 #ifdef FOO
00227 !     Control range and shape
00228 
00229       if (range(1) < shape (1) .or. shape (2) < range (2)) then
00230          print *, "range", range
00231          print *, "shape", shape
00232          call psmile_assert (__FILE__, __LINE__, &
00233                              "range must be within shape")
00234       endif
00235 
00236 !     Control control and range
00237 
00238       if (control(1) < range (1) .or. range (2) < control (2)) then
00239          print *, "control", control
00240          print *, "range  ", range
00241          call psmile_assert (__FILE__, __LINE__, &
00242                              "control must be within range")
00243       endif
00244 #endif /* FOO */
00245 #endif
00246 !
00247 ! temporaeren speicherbereich (control (2,1)-control(1,1)+1)*ndtry anlegen
00248 ! und do n =1,nref vor i-loop ziehen ?!?
00249 !
00250 ! ? statt immer alle nref zu kontrollieren, einmal location in der
00251 ! Zelle bestimmen ?
00252 !
00253 ! ? ein-dimensionaler Loop statt i,j loop
00254 !
00255 #ifdef DEBUG
00256       nfnd (:) = 0
00257 #endif /* DEBUG */
00258 
00259 #ifdef DEBUGX
00260       print *, 'begin ---'
00261         do k = control (1, 3), control (2,3)   
00262             do j = control(1,2), control (2,2) 
00263                do i = control(1,1), control (2,1) 
00264                print *, 'i,j,k,coord, loc, found', i,j,k, coords1 (i,j,k), &
00265                         loc(i,j,k), found (i,j,k)
00266                end do
00267             end do
00268         end do
00269 #endif
00270 !
00271 !===>
00272 !
00273         do k = control (1, 3), control (2,3)   
00274             do j = control(1,2), control (2,2) 
00275 !
00276 !-----------------------------------------------------------------------
00277 !             Look in the fine parts of the coarse grid cell
00278 !             loc(i)
00279 !
00280 !  Besser ein gatter/Scatter hier einbauen
00281 !-----------------------------------------------------------------------
00282 !
00283             ibegl = control (1, 1)
00284 !
00285 !===> ... Look for the next cells to be controlled
00286 !
00287             do while (ibegl <= control(2,1))
00288 !cdir vector
00289                   do i = ibegl, control(2,1)
00290                   if (found(i, j, k) == lev) exit
00291                   end do
00292 !
00293             if (i > control(2,1)) exit
00294             ibegl = i
00295 !
00296 !cdir vector
00297                do i = ibegl, control(2,1)
00298                if (found(i, j, k) /= lev) exit
00299                end do
00300 !
00301             iendl = i - 1
00302 !
00303 !===> ... Look in the real cell (ic)
00304 !
00305 #ifdef PRISM_ASSERTION
00306             do i = ibegl, iendl
00307             if (loc (i,j,k) < coords_shape(1) .or. &
00308                 loc (i,j,k) > coords_shape(2)) exit
00309             end do
00310 
00311             if (i <= iendl) then
00312                print *, "i,j,k, ic, coords_shape", &
00313                         i,j,k, loc (i,j,k), coords_shape
00314 print *, 'control', control
00315 print *, 'range  ', range
00316                call psmile_assert (__FILE__, __LINE__, &
00317                                    'wrong index')
00318             endif
00319 #endif
00320 !
00321             do jpart = 0, (iendl-ibegl)/5000
00322 !cdir vector loopcnt=5000
00323             do i = ibegl+jpart*5000, min (iendl, ibegl+(jpart+1)*5000 - 1)
00324 !
00325             ic = loc (i,j,k)
00326 !
00327             irange (1) = max (ic-1, grid_valid_shape (1))
00328             irange (2) = min (ic+1, grid_valid_shape (2))
00329 !
00330 !===> ... Compute distance to points
00331 !
00332                do n = irange(1), irange (2)
00333                dist (n-irange(1)+1) = abs (x_coords (n)-coords1(i,j,k))
00334                end do
00335 !
00336 !===> ... Look for the minimum distance
00337 !
00338             nref = irange (2) - irange (1) + 1
00339 #ifdef PRISM_ASSERTION
00340             if (nref <= 0) then
00341                call psmile_assert (__FILE__, __LINE__, &
00342                                    'nref <= 0')
00343             endif
00344 #endif
00345             nmin = MINLOC (dist(1:nref))
00346 #ifdef MINLOCFIX
00347             if (nmin(1).eq.0) nmin=1
00348 #endif /* MINLOCFIX */
00349 !
00350             loc (i,j,k) = irange (1) + nmin(1) - 1
00351 !
00352 #ifdef DEBUG2
00353    print *, "psmile_mg_method_1d_dble: i, j, k --- ", i,j,k
00354    print *, "psmile_mg_method_1d_dble loc   ", loc(i,j,k), nmin(1), dist(nmin(1))
00355    print *, "psmile_mg_method_1d_dble ic    ", ic
00356    print *, "psmile_mg_method_1d_dble irange", irange
00357    print *, "psmile_mg_method_1d_dble search", coords1(i,j,k)
00358    print *, "psmile_mg_method_1d_dble source", x_coords (irange(1):irange(2))
00359 #endif
00360 !
00361             if (dist(nmin(1)) .lt. tol) then
00362 !
00363 !===> ...... Point is located on the grid point of the method.
00364 !            Mark point by negative level number
00365 !
00366 !              found (i,j,k) = val_direct ! == lev == 1
00367 !
00368 #ifdef DEBUG
00369                nfnd (1) = nfnd (1) + 1
00370 #endif /* DEBUG */
00371             else
00372                found (i,j,k) = val_coupler
00373 #ifdef DEBUG
00374                nfnd (2) = nfnd (2) + 1
00375 #endif /* DEBUG */
00376 
00377                if (coords1(i,j,k) < chmin (loc(i,j,k)) .or. &
00378                    coords1(i,j,k) > chmax (loc(i,j,k))) then
00379 
00380                   irange (1) = max (ic-1, grid_valid_shape (1)-1)
00381 
00382                   do ii = irange (1), irange (2)
00383                   if (chmin(ii) <= coords1(i,j,k) .and. &
00384                       chmax(ii) >= coords1(i,j,k)) exit
00385                   end do
00386 !
00387                   if (ii <= irange(2)) then
00388 !
00389                      loc (i,j,k) = ii
00390                   else
00391 !                    Point as not found in the method grid;
00392 !                    mark point as 0
00393 #ifdef DEBUG
00394                      nfnd (3) = nfnd (3) + 1
00395                      nfnd (2) = nfnd (2) - 1
00396 #endif /* DEBUG */
00397 !Todo                found (i,j,k) = 0
00398  
00399                   endif
00400                endif
00401             endif
00402 !
00403             end do ! i
00404             end do ! jpart
00405 !
00406 !===> ... Is there an cyclic index
00407 !
00408             if (cyclic) then
00409 !
00410 !              This is an cyclic index;
00411 !              Move points, which are located in virtual cell 0
00412 !              into cell "grid_valid_shape(2)"
00413 !
00414 !cdir vector
00415                   do i = ibegl, iendl
00416                   if (loc (i,j,k) < grid_valid_shape (1)) then
00417                       loc (i,j,k) = grid_valid_shape (2)
00418                   endif
00419                   enddo
00420 #ifdef REMOVED
00421             else
00422 !
00423 !              Move points located in virtual cell 0 into
00424 !              cell 1; wird derzeit nur gemacht, um Probleme
00425 !              in anderen Routinen zu verhindern
00426 !              spaeter: diese indices in extra_search aufnehmen
00427 !              oder nearest neighbour
00428 !cdir vector
00429                   do i = ibegl, iendl
00430                   loc (i,j,k) = max (loc(i,j,k), grid_valid_shape (1))
00431                   enddo
00432 #endif
00433 !
00434             endif
00435 !
00436 !===> ... Start index of the next search
00437 !
00438             ibegl = iendl + 2
00439             end do ! while
00440 !
00441             end do ! j
00442          end do ! k
00443 !
00444 !===> All done
00445 !
00446 #ifdef DEBUGX
00447       print *, 'end ---'
00448         do k = control (1, 3), control (2,3)   
00449             do j = control(1,2), control (2,2) 
00450                do i = control(1,1), control (2,1) 
00451                print *, 'i,j,k,coord, loc, found', i,j,k, coords1 (i,j,k), &
00452                         loc(i,j,k), found (i,j,k)
00453                end do
00454             end do
00455         end do
00456 #endif
00457 
00458 #ifdef DEBUG
00459       print 9970, trim(ch_id), lev, nfnd
00460 9970  format (1x, a, ': psmile_mg_method_1d_dble: lev =', i3, &
00461               ': not :', i5, ', dir: ', i6, ', fnd :', i6, i5)
00462 #endif /* DEBUG */
00463 
00464 #ifdef VERBOSE
00465       print 9980, trim(ch_id), lev
00466 
00467       call psmile_flushstd
00468 #endif /* VERBOSE */
00469 !
00470 !  Formats:
00471 !
00472 
00473 9990 format (1x, a, ': psmile_mg_method_1d_dble: level', i3, &
00474                     ', control', 6i6)
00475 9980 format (1x, a, ': psmile_mg_method_1d_dble: eof, level', i3)
00476 
00477       end subroutine psmile_mg_method_1d_dble

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1