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

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1