psmile_mg_method_irreg2_real.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_irreg2_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_mg_method_irreg2_real (comp_info,    &
00012                       found, locations, search,              &
00013                       array, shape_2d, range_2d, control_2d, &
00014                              shape_1d, range_1d, control_1d, &
00015                       m_arrays, m_levdim,                    &
00016                       grid_id, method_id, tol, ierror)
00017 !
00018 ! !USES:
00019 !
00020       use PRISM_constants
00021       use psmile_grid, only : common_grid_range
00022       use PSMILe, dummy_interface => PSMILe_MG_method_irreg2_real
00023 
00024       Implicit none
00025 !
00026 ! !INPUT PARAMETERS:
00027 !
00028       Type (Enddef_comp), Intent (In) :: comp_info
00029 
00030 !     Info on the component in which the donor cells
00031 !     should be searched.
00032 
00033       Type (Enddef_search)            :: search
00034 
00035 !     Info's on coordinates to be searched
00036 
00037       Integer, Intent (In)            :: grid_id
00038 !
00039 !     Grid id
00040 
00041       Integer, Intent (In)            :: method_id
00042 !
00043 !     Method id
00044 
00045       Type (real_vector), Intent (In) :: array (ndim_3d, search%npart)
00046 !
00047       Integer,            Intent (In) :: shape_2d   (2, ndim_3d, search%npart)
00048       Integer,            Intent (In) :: range_2d   (2, ndim_3d, search%npart)
00049       Integer,            Intent (In) :: control_2d (2, ndim_3d, search%npart)
00050 !
00051       Integer,            Intent (In) :: shape_1d   (2, ndim_3d, search%npart)
00052       Integer,            Intent (In) :: range_1d   (2, ndim_3d, search%npart)
00053       Integer,            Intent (In) :: control_1d (2, ndim_3d, search%npart)
00054 
00055       Type (Enddef_mg_real)           :: m_arrays
00056       Integer,            Intent (In) :: m_levdim (ndim_3d)
00057 
00058       Real,   Intent (In)             :: tol
00059 !
00060 !     Absolute tolerance for search of "identical" points
00061 !     TOL >= 0.0
00062 !
00063 ! !INPUT/OUTPUT PARAMETERS:
00064 !
00065       Type (integer_vector)           :: found (search%npart, 2)
00066 
00067 !     found (*, 1) = Finest level number
00068 !     on which a grid cell was found in the irregular (corner) grid
00069 !     for point I.
00070 !
00071 !     found (*, 2) = Finest level number
00072 !     on which a grid cell was found in the regular (corner) grid
00073 !     for point I.
00074 !     Level number = -(nlev+1): Never found (input value)
00075 !
00076       Type (integer_vector)           :: locations (search%npart, 2)
00077 !
00078 !     locations (*, 1) = Indices of the grid cell (in first 2 dimensions)
00079 !     in which the point was found in the irregular (corner) grid.
00080 !
00081 !     locations (*, 2) = Indices of the grid cell (3rd direction)
00082 !     in which the point was found in the regular (corner) grid.
00083 !
00084 ! !OUTPUT PARAMETERS:
00085 !
00086       Integer, Intent (Out)           :: ierror
00087 
00088 !     Returns the error code of PSMILe_MG_method_irreg2_real;
00089 !             ierror = 0 : No error
00090 !             ierror > 0 : Severe error
00091 !
00092 ! !LOCAL VARIABLES
00093 !
00094 ! indl = Index of LonLat values in "found, locations"
00095 ! indz = Index of Vert   values in "found, locations"
00096 ! nc_reg = Number of corners for regular directions
00097 !
00098       Integer, Parameter            :: indl = 1
00099       Integer, Parameter            :: indz = 2
00100       Integer, Parameter            :: nc_reg = 2
00101 !
00102       Integer, Parameter            :: lev = 1
00103 !
00104 !     ... for temporay storage
00105 !
00106       Integer                       :: i, ipart
00107 !
00108 !     ... for grid and methods
00109 !
00110       Type(Grid),          Pointer  :: grid_info
00111       Type (Coords_Block), Pointer  :: coords_pointer
00112       Type (Corner_Block), Pointer  :: corner_pointer
00113       Integer,             Pointer  :: grid_valid_shape (:, :)
00114 !
00115 !     ... for virtual cells
00116 !
00117       Integer,          Allocatable :: bmaski (:, :)
00118       Integer,          Allocatable :: bmaskj (:, :)
00119 !
00120       Integer                       :: c_levdim (2, ndim_2d)
00121       Type (Enddef_mg_real), Pointer   :: c_arrays
00122 !
00123 !     ... for cyclic cells
00124 !
00125       Real, Save                    :: period (ndim_2d)
00126 !
00127 !     ... for error parameters
00128 !
00129       Integer, Parameter            :: nerrp = 2
00130       Integer                       :: ierrp (nerrp)
00131 !
00132 ! !DESCRIPTION:
00133 !
00134 ! Subroutine "PSMILe_MG_method_irreg2_real" searches the donor cells
00135 ! for the subgrid sent by the sending process, if
00136 !   the first 2 grid coordinates are irregular and
00137 !   the last    grid coordinate  is    regular.
00138 !
00139 ! Subroutine "PSMILe_MG_method_irreg2_real" is designed for grids of
00140 ! type "PRISM_Irrlonlat_regvrt". We assume that this routine is only called
00141 ! for methods that are used for coupling.
00142 !
00143 !
00144 ! !REVISION HISTORY:
00145 !   Date      Programmer   Description
00146 ! ----------  ----------   -----------
00147 ! 05.05.09    H. Ritzdorf  created
00148 !
00149 !EOP
00150 !----------------------------------------------------------------------
00151 !
00152 ! $Id: psmile_mg_method_irreg2_real.F90 2833 2010-12-21 14:14:42Z coquart $
00153 ! $Author: coquart $
00154 !
00155    Character(len=len_cvs_string), save :: mycvs = 
00156        '$Id: psmile_mg_method_irreg2_real.F90 2833 2010-12-21 14:14:42Z coquart $'
00157 !
00158 !----------------------------------------------------------------------
00159 !
00160 !  Initialization
00161 !
00162 #ifdef VERBOSE
00163       print 9990, trim(ch_id), Grids(grid_id)%comp_id
00164 
00165       call psmile_flushstd
00166 #endif /* VERBOSE */
00167 
00168 #ifdef PRISM_ASSERTION
00169       if (method_id > Number_of_Methods_allocated .or. &
00170           method_id < 1) then
00171          print *, trim(ch_id), "PSMILe_MG_method_irreg2_real: method_id =", &
00172                   method_id
00173          call psmile_assert (__FILE__, __LINE__, &
00174               "Invalid Method id")
00175       endif
00176 #endif
00177 
00178 !     Initialization
00179 
00180       period (1:ndim_2d) = common_grid_range(2,1:ndim_2d) - &
00181                            common_grid_range(1,1:ndim_2d)
00182 
00183       grid_info      => Grids(grid_id)
00184       corner_pointer => grid_info%corner_pointer
00185       coords_pointer => Methods(method_id)%coords_pointer
00186 !
00187       grid_valid_shape  => grid_info%grid_shape
00188 !
00189       Allocate (bmaski (grid_valid_shape(1,1)-1:grid_valid_shape(2,1)+1, 2), &
00190                 bmaskj (grid_valid_shape(1,2)-1:grid_valid_shape(2,2)+1, 2), &
00191                STAT = ierror)
00192       if ( ierror > 0 ) then
00193          ierrp (1) = ierror
00194          ierrp (2) = 2*(grid_valid_shape(2,1)-grid_valid_shape(1,1)+1+2)*2
00195          ierror = PRISM_Error_Alloc
00196          call psmile_error ( ierror, 'bmask{i,j}', &
00197                              ierrp, 2, __FILE__, __LINE__ )
00198          return
00199       endif
00200 
00201 !     ... for the virtual cells of the 2-dimensional block
00202 
00203       c_arrays => grid_info%mg_infos(lev)%real_arrays
00204 !
00205       c_levdim (1, 1:ndim_2d) = grid_info%ijk0 (1:ndim_2d)
00206       c_levdim (2, 1:ndim_2d) = c_levdim (1, 1:ndim_2d) + &
00207          grid_info%mg_infos(lev)%levdim(1:ndim_2d)
00208 
00209       ! Rene: cyclic and period are not need in
00210       ! PSMILe_bbcells_virt_2d_real and can be erased.
00211 
00212       call psmile_bbcells_virt_2d_real ( method_id,                  &
00213                 coords_pointer%coords_real(1)%vector,                &
00214                 coords_pointer%coords_real(2)%vector,                &
00215                 coords_pointer%coords_shape, grid_valid_shape,       &
00216                 corner_pointer%corners_real(1)%vector,               &
00217                 corner_pointer%corners_real(2)%vector,               &
00218                 corner_pointer%corner_shape(:,1:ndim_2d),            &
00219                 corner_pointer%nbr_corners/nc_reg,                   &
00220                 c_arrays%chmin(1)%vector, c_arrays%chmin(2)%vector,  &
00221                 c_arrays%chmax(1)%vector, c_arrays%chmax(2)%vector,  &
00222                 c_levdim,                                            &
00223                 m_arrays%chmin(1)%vector, m_arrays%chmin(2)%vector,  &
00224                 m_arrays%chmax(1)%vector, m_arrays%chmax(2)%vector,  &
00225                 m_arrays%midp(1)%vector,  m_arrays%midp(2)%vector,   &
00226                 m_levdim, period (1), bmaski, bmaskj, ierror)
00227       if (ierror > 0) return
00228 !
00229 !     Compute bounding boxes for the cells
00230 !
00231 !     ... within the 2-dimensional block
00232 !
00233       do i = 1, ndim_2d
00234       call psmile_bbcells_2d_real (                                 &
00235                   coords_pointer%coords_real(i)%vector,                &
00236                   coords_pointer%coords_shape, grid_valid_shape,       &
00237                   corner_pointer%corner_shape(:,1:ndim_2d),            &
00238                   m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector,  &
00239                   m_arrays%midp(i)%vector, m_levdim,                   &
00240                   period (i), ierror)
00241       if (ierror > 0) return
00242       end do
00243 !
00244 !     ... for the z-direction
00245 !
00246       i = ndim_3d
00247       call psmile_bbcells_1d_real (                                     &
00248                 coords_pointer%coords_real(i)%vector,                   &
00249                 coords_pointer%coords_shape(1:2, i),                    &
00250                 grid_valid_shape  (1:2, i),                             &
00251                 corner_pointer%corners_real(i)%vector,                  &
00252                 corner_pointer%corner_shape(1:2,i), nc_reg,             &
00253                 m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector,     &
00254                 m_levdim(i), ierror)
00255       if (ierror > 0) return
00256 !
00257 !     ... correction of pole cells
00258 !
00259       if (associated(corner_pointer%pole_array)) then
00260          call psmile_bbcells_pole_real (                              &
00261                   coords_pointer%coords_shape,                        &
00262                   coords_pointer%coords_real(1)%vector,               &
00263                   coords_pointer%coords_real(2)%vector,               &
00264                   corner_pointer%corner_shape,                        &
00265                   grid_valid_shape,                                   &
00266                   m_arrays%chmin(1)%vector, m_arrays%chmax(1)%vector, &
00267                   m_arrays%chmin(2)%vector, m_arrays%chmax(2)%vector, &
00268                   m_arrays%midp(1)%vector, m_arrays%midp(2)%vector,   &
00269                   corner_pointer%pole_array, period(1), ierror)
00270          if (ierror > 0) return
00271       endif
00272 
00273 !
00274 !===> Get locations in first 2 irregular directions
00275 !
00276          do ipart = 1, search%npart
00277 
00278 #ifdef DEBUG
00279            PRINT*, TRIM(ch_id), ' intersection :',ipart
00280            call psmile_flushstd
00281 #endif
00282 
00283          call psmile_mg_method_2d_real ( comp_info, grid_info%nlev,        &
00284                    found(ipart,indl)%vector, locations(ipart,indl)%vector, &
00285                    range_2d (:, :, ipart),                                 &
00286                    array(1,ipart)%vector, array(2,ipart)%vector,           &
00287                    shape_2d (:, :, ipart), control_2d (:, :, ipart),       &
00288                    coords_pointer%coords_real(1)%vector,                   &
00289                    coords_pointer%coords_real(2)%vector,                   &
00290                    coords_pointer%coords_shape, grid_valid_shape,          &
00291                    grid_info%cyclic, period,                               &
00292                    m_arrays%chmin(1)%vector, m_arrays%chmin(2)%vector,     &
00293                    m_arrays%chmax(1)%vector, m_arrays%chmax(2)%vector,     &
00294                    tol, ierror)
00295          if (ierror > 0) return
00296 !
00297 !===> Get locations in 3rd regular direction
00298 !
00299          call psmile_mg_method_1d_real ( comp_info, grid_info%nlev,        &
00300                    found(ipart,indz)%vector, locations(ipart,indz)%vector, &
00301                    range_1d (:, :, ipart),                                 &
00302                    search%search_real(ndim_3d, ipart)%vector,              &
00303                    shape_1d (:, :, ipart), control_1d (:, :, ipart),       &
00304                    method_id,                                              &
00305                    coords_pointer%coords_real(ndim_3d)%vector,             &
00306                    coords_pointer%coords_shape(:, ndim_3d),                &
00307                    grid_valid_shape  (:, ndim_3d),                         &
00308                    grid_info%cyclic(ndim_3d),                              &
00309                    m_arrays%chmin(ndim_3d)%vector,                         &
00310                    m_arrays%chmax(ndim_3d)%vector, tol, ierror)
00311               if (ierror > 0) return
00312 !
00313          end do ! ipart
00314 !
00315 !===> Deallocate temporary arrays
00316 !
00317       Deallocate (bmaski, bmaskj)
00318 
00319 #ifdef VERBOSE
00320       print 9980, trim(ch_id), grid_info%comp_id, ierror
00321 
00322       call psmile_flushstd
00323 #endif /* VERBOSE */
00324 !
00325 !  Formats:
00326 !
00327 9990 format (1x, a, ': psmile_mg_method_irreg2_real: comp_id =', i3)
00328 9980 format (1x, a, ': psmile_mg_method_irreg2_real: comp_id =', i3, &
00329                     '; eof ierror =', i4)
00330 
00331       end subroutine PSMILe_MG_method_irreg2_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1