psmile_locations_irreg2.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_Locations_irreg2
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_locations_irreg2 (                   &
00012                       found, loc, range, control,            &
00013                       search, method_id, msk_required,       &
00014                       virtual_cell, virtual_cell_required,   &
00015                       dir_index, cpl_index, len_cpl, ierror)
00016 !
00017 ! !USES:
00018 !
00019       use PRISM_constants
00020 !
00021       use PSMILe, dummy_interface => PSMILe_Locations_irreg2
00022 
00023       implicit none
00024 !
00025 ! !INPUT PARAMETERS:
00026 !
00027       Type (Enddef_search), Intent (InOut) :: search
00028 
00029 !     Info's on coordinates to be searched
00030 !
00031       Integer, Intent (In)            :: range (2, ndim_3d, search%npart)
00032 
00033 !     Dimension of loc and found
00034 !     ??? muss range nicht shape sein ?
00035 !
00036       Type (integer_vector)           :: found (search%npart, 2)
00037 
00038 !     Was a grid cell found for for point (i,j,k) in range
00039 !     (2, ndim_3d, ipart) for ipart-th partition ?
00040 !     Only found (1:2, *) is significant.
00041 !     found(*,i)%vector(n) = val_direct : Point was found and
00042 !                                         has to be sent directly to target 
00043 !                                         process.
00044 !     found(*,i)%vector(n) = val_coupler: Point was found and
00045 !                                         has to be sent to coupler process
00046 !     Otherwise                           Point was not found.
00047 
00048       Type (integer_vector)           :: loc (search%npart, 2)
00049 !
00050 !     Indices of the grid cell in which the point was found.
00051 !     Assumed input value loc(:)%vector(:, len) = 0
00052 !     Only loc (1:2, *) is significant.
00053 
00054       Type (integer_vector)           :: virtual_cell (search%npart)
00055 
00056 !     Code for virtual cells if point was found in virtual cell
00057 
00058       Integer, Intent (In)            :: control (2, ndim_3d, search%npart)
00059 
00060 !     Index range found
00061 !
00062       Integer, Intent (In)            :: method_id
00063 
00064 !     Method id of donor cells (on this process)
00065 
00066       Logical, Intent (In)            :: msk_required
00067 
00068 !     Switch to determine whether an additional src location mask is required.
00069 !
00070       Logical, Intent (In)            :: virtual_cell_required
00071 
00072 !     Storage of virtual cell info required ?
00073 
00074 ! !OUTPUT PARAMETERS:
00075 !
00076       Integer, Intent (Out)           :: dir_index
00077 
00078 !     Send index if data has to be directly transferred to the
00079 !     destination process.
00080 !     PRISM_undefined, if no data has to be sent.
00081 !
00082       Integer, Intent (Out)           :: cpl_index
00083 
00084 !     Send index if data has to be transferred to the coupler
00085 !     PRISM_undefined, if no data has to be sent.
00086 
00087       Integer, Intent (Out)           :: len_cpl (search%npart)
00088 
00089 !     Number of points to be sent to the coupler for each partition.
00090 !
00091       Integer, Intent (Out)           :: ierror
00092 
00093 !     Returns the error code of PSMILe_Locations_irreg2;
00094 !             ierror = 0 : No error
00095 !             ierror > 0 : Severe error
00096 !
00097 ! !LOCAL PARAMETERS
00098 !
00099 ! indl = Index of LonLat values in "found, loc, ibeg, iend"
00100 ! indz = Index of Vert   values in "found, loc, ibeg, iend"
00101 !
00102       Integer, Parameter              :: val_direct  =  1
00103       Integer, Parameter              :: val_coupler = -1
00104       Integer, Parameter              :: val_both    =  0
00105 !
00106       Integer, Parameter              :: send_coupler_index = 1
00107       Integer, Parameter              :: send_direct_index  = 2
00108 !
00109       Integer, Parameter              :: indl = 1
00110       Integer, Parameter              :: indz = 2
00111 ! !LOCAL VARIABLES
00112 !
00113 !     ... Method pointer
00114 !
00115       Type (Method), Pointer          :: mp
00116 !
00117 !     ... for locations to be stored and transferred
00118 !
00119       Real                            :: ratio
00120 !
00121       Integer                         :: ncpl_tot, ndir_tot
00122       Integer                         :: ncpl  (search%npart)
00123       Integer                         :: ndir  (search%npart)
00124       Integer                         :: ncplz (search%npart)
00125       Integer                         :: ndirz (search%npart)
00126       Integer                         :: index, n, val_cpl
00127 !
00128 !     ... for partitions
00129 !
00130       Integer                         :: ipart, nadd, nprev
00131       Integer                         :: ibeg (search%npart, indz)
00132       Integer                         :: iend (search%npart, indz)
00133       Integer                         :: nprevi
00134 !
00135 !     ... for error handling
00136 !
00137 !     Integer, parameter              :: nerrp = 3
00138 !     Integer                         :: ierrp (nerrp)
00139 !
00140 ! !DESCRIPTION:
00141 !
00142 ! Subroutine "PSMILe_Locations_irreg2" stores the data on
00143 ! locations found for the method (grid) and the subgrid coords.
00144 !
00145 ! Subroutine "PSMILe_Locations_irreg2" is designed for grids of
00146 ! of type "PRISM_Irrlonlat_regvrt".
00147 !
00148 !
00149 ! !REVISION HISTORY:
00150 !
00151 !   Date      Programmer   Description
00152 ! ----------  ----------   -----------
00153 ! 03.07.04    H. Ritzdorf  created
00154 !
00155 !EOP
00156 !----------------------------------------------------------------------
00157 !
00158 !  $Id: psmile_locations_irreg2.F90 2938 2011-02-03 10:23:22Z hanke $
00159 !  $Author: hanke $
00160 !
00161    Character(len=len_cvs_string), save :: mycvs = 
00162        '$Id: psmile_locations_irreg2.F90 2938 2011-02-03 10:23:22Z hanke $'
00163 !
00164 !----------------------------------------------------------------------
00165 !
00166 !  Initialization
00167 !
00168 #ifdef VERBOSE
00169       print 9990, trim(ch_id), method_id, search%msg_intersections%first_tgt_method_id, &
00170                                search%sender
00171 
00172       call psmile_flushstd
00173 #endif /* VERBOSE */
00174 !
00175       ierror = 0
00176 !
00177       cpl_index = PRISM_undefined
00178       dir_index = PRISM_undefined
00179 !
00180       len_cpl = 0
00181 !
00182       mp => Methods(method_id)
00183 !
00184 #ifdef PRISM_ASSERTION
00185       if (search%grid_type == PRISM_Irrlonlatvrt) then
00186          print *, "search%grid_type = ", search%grid_type
00187          call psmile_assert (__FILE__, __LINE__, &
00188                              "This routine is not designed for this grid type")
00189       endif
00190 !
00191       if ( mp%status == PSMILe_Status_free .or. &
00192            mp%status == PSMILe_Status_undefined ) then
00193          call psmile_assert (__FILE__, __LINE__, "Incorrect method")
00194       endif
00195 
00196          do ipart = 1, search%npart
00197          if (control (1,1,ipart) <   range (1,1,ipart) .or. &
00198                range (2,1,ipart) < control (2,1,ipart) .or. &
00199              control (1,2,ipart) <   range (1,2,ipart) .or. &
00200                range (2,2,ipart) < control (2,2,ipart) .or. &
00201              control (1,3,ipart) <   range (1,3,ipart) .or. &
00202                range (2,3,ipart) < control (2,3,ipart)) then
00203             print *, ipart, control (:, :, ipart), range (:, :, ipart)
00204             call psmile_assert (__FILE__, __LINE__, "control is out of range")
00205          endif
00206          end do ! ipart
00207 #endif
00208 !
00209 !----------------------------------------------------------------------------
00210 !     Determine maximal number of points to be transferred to the
00211 !     coupler or directly to the destination process
00212 !----------------------------------------------------------------------------
00213 !
00214       ndir_tot = 0
00215       ncpl_tot = 0
00216 !
00217       ndir(:) = 0
00218       ncpl(:) = 0
00219 !
00220       ndirz(:) = 0
00221       ncplz(:) = 0
00222 !
00223       ibeg (:, :) = 1
00224 !
00225          do ipart = 1, search%npart
00226 !
00227 !===> Get number of locations to be sent directly to another process
00228 !     ndir = Number of locations which can be sent directly
00229 !     ncpl = Number of locations sent to the coupler
00230 !
00231 !     ? ibeg_dir, iend_dir, ibeg_cpl, iend_cpl berechnen ?
00232 !
00233 !     ibeg = (control(1,2)-range(1,2))*(range(2,1)-range(1,1)+1) &
00234 !          + (control(1,1)-range(1,1)) + 1
00235 !
00236 !     iend = (control(2,2)-range(1,2))*(range(2,1)-range(1,1)+1) &
00237 !          + (control(2,1)-range(1,1)) + 1
00238 !
00239 !     ... Get number of locations in 2d fields "found" and "loc"
00240 !
00241          iend (ipart, indl) = (range (2,1, ipart) - range(1,1, ipart) + 1) * &
00242                               (range (2,2, ipart) - range(1,2, ipart) + 1)
00243          iend (ipart, indz) =  range (2,3, ipart) - range(1,3, ipart) + 1
00244 !
00245 !cdir vector
00246             do n = ibeg(ipart, indl), iend (ipart, indl)
00247             if (found(ipart,indl)%vector(n) == val_coupler) then
00248                ncpl(ipart) = ncpl(ipart) + 1
00249             else if (found(ipart,indl)%vector(n) == val_direct) then
00250                ndir(ipart) = ndir(ipart) + 1
00251             endif
00252             end do
00253 !
00254 !     ... Get number of locations in 1d fields "foundz" and "locz"
00255 !
00256 !cdir vector
00257             do n = ibeg(ipart, indz), iend (ipart, indz)
00258             if (found(ipart,indz)%vector(n) == val_coupler) then
00259                ncplz(ipart) = ncplz(ipart) + 1
00260             else if (found(ipart,indz)%vector(n) == val_direct) then
00261                ndirz(ipart) = ndirz(ipart) + 1
00262             endif
00263             end do
00264 
00265          end do ! ipart
00266 
00267 #ifdef DEBUG
00268       print *, trim(ch_id), ': psmile_locations_irreg2:'
00269       print *, "range:  ", range
00270       print *, "control:", control
00271       print *, "ibeg:",    ibeg
00272       print *, "iend:",    iend
00273       print *, "ncpl   :", ncpl, ncplz
00274       print *, "ndir   :", ndir, ndirz
00275 
00276       call psmile_flushstd
00277 #endif
00278 !
00279 !     ... Get total number of locations
00280 
00281 ! Was ist mit search%grid_type ?
00282 
00283          do ipart = 1, search%npart
00284          n = ndir(ipart) * ndirz(ipart)
00285          ndir_tot = ndir_tot + n
00286          ncpl_tot = ncpl_tot + ((ncpl (ipart)+ndir (ipart)) * &
00287                                 (ncplz(ipart)+ndirz(ipart)) - n)
00288          end do ! ipart
00289 !
00290       if (ncpl_tot + ndir_tot > 0) then
00291          ratio = real(ndir_tot) / (ncpl_tot+ndir_tot)
00292       else
00293          ratio = 0.0
00294       endif
00295 !
00296 #ifdef ONLY_FOR_TESTING
00297       print *, '######## psmile_locations_irreg2: ratio set to 0.0'
00298       ratio = 0.01
00299 #endif
00300 #ifdef DEBUG
00301       print *, 'ncpl_tot, ndir_tot, ratio:', ncpl_tot, ndir_tot, ratio 
00302 #endif
00303 !
00304 !===> Store locations
00305 !     ??? Right now data a stored row by row.
00306 !         Clustering of Data, if possible !!!
00307 !     ??? ratio is  weak criteria;
00308 !
00309       if (max(maxval (ndir), maxval (ndirz)) > 0 .and. ratio < 0.05) then
00310 !
00311 !       ... The number of points which can be sent directly is not
00312 !           "large" enough
00313 !       ... Send all points to the coupler
00314 !
00315         ncpl = ncpl + ndir
00316         ndir = 0
00317 !
00318         ncplz = ncplz + ndirz
00319         ndirz = 0
00320 !
00321         ncpl_tot = ncpl_tot + ndir_tot
00322         ndir_tot = 0
00323 !
00324         val_cpl = val_both
00325       else
00326         val_cpl = val_coupler
00327 !       ... make vertical points for direct send also available for coupling
00328         ncplz = ncplz + ndirz
00329       endif
00330 !
00331 !----------------------------------------------------------------------------
00332 !     Generate send areas for data send to a coupler process
00333 !----------------------------------------------------------------------------
00334 !
00335       if (ncpl_tot > 0) then
00336         call psmile_get_info_index (method_id, send_coupler_index, index, ierror)
00337         if (ierror > 0) return
00338 
00339         cpl_index = index
00340 !
00341 !       ... generate send areas for data exchange with coupler
00342 !
00343 !       mp%send_infos_coupler(index)%nloc = ncpl_tot
00344         mp%send_infos_coupler(index)%nvec = 2
00345         mp%send_infos_coupler(index)%nparts = search%npart
00346 !
00347         mp%send_infos_coupler(index)%nrecv    = 0
00348         mp%send_infos_coupler(index)%num2recv = 0
00349 !
00350         mp%send_infos_coupler(index)%remote_method_id = search%msg_intersections%first_tgt_method_id
00351 !
00352 !       ... Rank of coupler in communicator "comm_coupler"
00353 !
00354         mp%send_infos_coupler(index)%dest = 0
00355 !
00356 !       ... generate send areas for data exchange with coupler
00357 !
00358         call psmile_locations_alloc (mp%send_infos_coupler(index), ierror)
00359         if (ierror > 0) return
00360 
00361         nprev  = 0
00362         nprevi = 0
00363 
00364          if (virtual_cell_required) then
00365             Allocate (mp%send_infos_coupler(index)%virtual(1,search%npart), &
00366                      STAT = ierror)
00367                if ( ierror > 0 ) then
00368 
00369                   call psmile_error ( PRISM_Error_Alloc, 'mp%send_infos_coupler(index)%virtual', &
00370                         (/ierror, search%npart/), 2, __FILE__, __LINE__ )
00371                   ierror = PRISM_Error_Alloc
00372                   return
00373                endif
00374                ! Derefference pointer explicitly as some compilers
00375                ! set a reference by default
00376                do ipart = 1, search%npart
00377                   Nullify(mp%send_infos_coupler(index)%virtual(1,ipart)%vector)
00378                enddo
00379          endif
00380 !
00381            do ipart = 1, search%npart
00382 !
00383            call psmile_store_dest_locs_21d (                &
00384               found(ipart,indl)%vector, range(1:2, 1:ndim_3d, ipart), &
00385               control (1:2, 1:ndim_3d, ipart), found(ipart,indz)%vector, &
00386               mp%send_infos_coupler(index), ncpl_tot, &
00387               val_cpl, nprev, nadd, ierror)
00388            if (ierror > 0) return
00389 !
00390            nprev = nprev + nadd
00391            len_cpl (ipart) = nadd
00392 
00393            call psmile_store_source_locs_2d (                          &
00394                      found(ipart,indl)%vector, loc(ipart,indl)%vector, &
00395                      ibeg (ipart,indl),       iend(ipart,indl),        &
00396                      mp%send_infos_coupler(index), ncpl(ipart),        &
00397                      val_cpl, indl, ipart, nprevi, nadd, ierror)
00398            if (ierror > 0) return
00399 !
00400 !          ... val_both: TODO: make vertical points for direct send
00401 !                        also available for coupling
00402 !
00403            call psmile_store_source_locs_1d (             &
00404                      found(ipart,indz)%vector, loc(ipart,indz)%vector, &
00405                      ibeg (ipart,indz),       iend(ipart,indz),  &
00406                      mp%send_infos_coupler(index), ncplz(ipart), &
00407                      val_both, indz, ipart, nprevi, nadd, ierror)
00408            if (ierror > 0) return
00409 !
00410            if ( msk_required ) then
00411               call psmile_store_mask_locs_21d (          &
00412                         range   (1:2, 1:ndim_2d, ipart), &
00413                         control (1:2, 1:ndim_2d, ipart), &
00414                         found(ipart,indl)%vector,        &
00415                         range   (1:2,   ndim_3d, ipart), &
00416                         control (1:2,   ndim_3d, ipart), &
00417                         found(ipart,indz)%vector,        &
00418                         mp%send_infos_coupler(index),    &
00419                         ipart, ncpl(ipart), ncplz(ipart), ierror )
00420               if (ierror > 0) return
00421            endif
00422 !
00423            if ( virtual_cell_required ) then
00424               call psmile_store_source_virt_3d (                     &
00425                         found(ipart,indl)%vector,                    &
00426                         virtual_cell(ipart)%vector,                  &
00427                         ibeg (ipart,indl), iend(ipart,indl),         &
00428                         mp%send_infos_coupler(index), ncpl(ipart),   &
00429                         val_cpl, indl, ipart, nprevi, ierror)
00430               if (ierror > 0) return
00431            endif
00432 
00433            enddo
00434 !
00435         mp%send_infos_coupler(index)%nloc = nprev
00436       endif
00437 !
00438 !----------------------------------------------------------------------------
00439 !     Generate send areas for data exchange with application process
00440 !----------------------------------------------------------------------------
00441 !
00442       if (ndir_tot > 0) then
00443 !
00444         call psmile_get_info_index (method_id, send_direct_index, index, ierror)
00445         if (ierror > 0) return
00446 !
00447         dir_index = index
00448 !
00449         mp%send_infos_direct(index)%dest = search%sender
00450 !       mp%send_infos_direct(index)%nloc = ndir_tot
00451 !
00452         mp%send_infos_direct(index)%nvec      = 2
00453         mp%send_infos_direct(index)%nparts    = search%npart
00454 !
00455         mp%send_infos_direct(index)%nrecv    = 0
00456         mp%send_infos_direct(index)%num2recv = 0
00457 !
00458         mp%send_infos_direct(index)%remote_method_id = search%msg_intersections%first_tgt_method_id
00459 !
00460 !       ... generate send areas for direct data exchange
00461 !
00462         call psmile_locations_alloc (mp%send_infos_direct(index), ierror)
00463         if (ierror > 0) return
00464 !
00465         nprev  = 0
00466         nprevi = 0
00467 !
00468            do ipart = 1, search%npart
00469 
00470            call psmile_store_dest_locs_21d (                &
00471               found(ipart,indl)%vector, range(1:2, 1:ndim_3d, ipart), &
00472                control (1:2, 1:ndim_3d, ipart), found(ipart,indz)%vector, &
00473               mp%send_infos_direct(index), ndir_tot,        &
00474               val_direct, nprev, nadd, ierror)
00475            if (ierror > 0) return
00476 !
00477            nprev = nprev + nadd
00478 
00479            call psmile_store_source_locs_2d (                   &
00480               found(ipart,indl)%vector, loc(ipart,indl)%vector, &
00481               ibeg (ipart,indl),       iend(ipart,indl),        &
00482               mp%send_infos_direct(index), ndir(ipart),         &
00483               val_direct, indl, ipart, nprevi, nadd, ierror)
00484            if (ierror > 0) return
00485 !
00486            call psmile_store_source_locs_1d (                   &
00487               found(ipart,indz)%vector, loc(ipart,indz)%vector, &
00488               ibeg (ipart,indz),       iend(ipart,indz),        &
00489               mp%send_infos_direct(index), ndirz(ipart),        &
00490               val_direct, indz, ipart, nprevi, nadd, ierror)
00491            if (ierror > 0) return
00492 !
00493            end do ! ipart
00494 !
00495         mp%send_infos_direct(index)%nloc = nprev
00496       endif
00497 !
00498 !===> All done
00499 !
00500 #ifdef VERBOSE
00501       print 9980, trim(ch_id), ierror, dir_index, cpl_index
00502 
00503       call psmile_flushstd
00504 #endif /* VERBOSE */
00505 !
00506 !  Formats:
00507 !
00508 
00509 #ifdef VERBOSE
00510 
00511 9990 format (1x, a, ': psmile_locations_irreg2: method_id', i3, &
00512                     ' to ', i3, '(', i2, ')')
00513 9980 format (1x, a, ': psmile_locations_irreg2: eof ierror =', i3, &
00514                     ', dir_index =', i10, ', cpl_index =', i10)
00515 
00516 #endif /* VERBOSE */
00517 
00518       end subroutine PSMILe_Locations_irreg2

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1