psmile_locations_3d_mask.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_mask
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_locations_3d_mask (search, inter, shift, method_id, &
00012                                            dir_index, ierror)
00013 !
00014 ! !USES:
00015 !
00016       use PRISM_constants
00017 !
00018       use PSMILe, dummy_interface => PSMILe_Locations_3d_mask
00019 
00020       implicit none
00021 !
00022 ! !INPUT PARAMETERS:
00023 
00024       Type (Enddef_search), Intent (InOut) :: search
00025 
00026 !     Info's on coordinates to be searched
00027 !
00028       Integer, Intent (In)            :: inter (2, ndim_3d, search%npart)
00029 
00030 !     Common intersections found
00031 
00032       Integer, Intent (In)            :: shift (ndim_3d)
00033 
00034 !     Shift to transfer local indices into local indices
00035 !     of target grid with
00036 !     target_grid_index (:, i) = local index (:, i) - shift (:)
00037 !
00038       Integer, Intent (In)            :: method_id
00039 
00040 !     Method id of donor cells (on this process)
00041 
00042 ! !OUTPUT PARAMETERS:
00043 !
00044       Integer, Intent (Out)           :: dir_index
00045 
00046 !     Send index if data has to be directly transferred to the
00047 !     destination process.
00048 !     PRISM_undefined, if no data has to be sent.
00049 !
00050       Integer, Intent (Out)           :: ierror
00051 
00052 !     Returns the error code of PSMILe_Locations_3d_mask;
00053 !             ierror = 0 : No error
00054 !             ierror > 0 : Severe error
00055 !
00056 ! !LOCAL PARAMETERS
00057 !
00058       Integer, Parameter              :: send_direct_index  = 2
00059 !
00060 ! !LOCAL VARIABLES
00061 !
00062 !     ... Method pointer
00063 !
00064       Type (Method), Pointer          :: mp
00065 !
00066 !     ... for locations to be stored and transferred
00067 !
00068       Integer                         :: n, nadd, ndir, nprev
00069       Integer                         :: index
00070 !
00071 !     ... for partitions
00072 !
00073       Integer                         :: ipart
00074 !
00075 !     ... for error handling
00076 !
00077 !     Integer, parameter              :: nerrp = 3
00078 !     Integer                         :: ierrp (nerrp)
00079 !
00080 ! !DESCRIPTION:
00081 !
00082 ! Subroutine "PSMILe_Locations_3d_mask" stores the data on locations to
00083 ! be sent from the source grid to the target grid.
00084 ! Subroutine "PSMILe_Locations_3d_mask" is designed for grids of type
00085 ! PRISM_Gridless where the locations are only defined by the mask entries
00086 ! and all data is sent directly to the target process.
00087 !
00088 ! TODO: request into send_info and free dstijk
00089 ! TODO: elimination of multiple points in dstijk and srclocs
00090 !
00091 ! !REVISION HISTORY:
00092 !
00093 !   Date      Programmer   Description
00094 ! ----------  ----------   -----------
00095 ! 11.02.05    H. Ritzdorf  created
00096 !
00097 !EOP
00098 !----------------------------------------------------------------------
00099 !
00100 !  $Id: psmile_locations_3d_mask.F90 2788 2010-11-30 14:34:07Z hanke $
00101 !  $Author: hanke $
00102 !
00103    Character(len=len_cvs_string), save :: mycvs = 
00104        '$Id: psmile_locations_3d_mask.F90 2788 2010-11-30 14:34:07Z hanke $'
00105 !
00106 !----------------------------------------------------------------------
00107 !
00108 !  Initialization
00109 !
00110 #ifdef VERBOSE
00111       print 9990, trim(ch_id), method_id, search%msg_intersections%first_tgt_method_id, &
00112                                search%sender
00113 
00114       call psmile_flushstd
00115 #endif /* VERBOSE */
00116 !
00117       ierror = 0
00118 !
00119       dir_index = PRISM_undefined
00120 !
00121       mp => Methods(method_id)
00122 !
00123 #ifdef PRISM_ASSERTION
00124       if ( mp%status == PSMILe_Status_free .or. &
00125            mp%status == PSMILe_Status_undefined ) then
00126          call psmile_assert (__FILE__, __LINE__, "Incorrect method")
00127       endif
00128 
00129          do ipart = 1, search%npart
00130          if (inter (1,1,ipart) < search%range (1,1,ipart)   .or. &
00131                search%range (2,1,ipart) < inter (2,1,ipart) .or. &
00132              inter (1,2,ipart) < search%range (1,2,ipart)   .or. &
00133                search%range (2,2,ipart) < inter (2,2,ipart) .or. &
00134              inter (1,3,ipart) < search%range (1,3,ipart)   .or. &
00135                search%range (2,3,ipart) < inter (2,3,ipart)) then
00136              print *, ipart, inter (:,:,ipart), search%range (:,:,ipart)
00137              call psmile_assert (__FILE__, __LINE__, "inter is out of range")
00138          endif
00139          end do ! ipart
00140 #endif
00141 !
00142 !----------------------------------------------------------------------------
00143 !     Determine maximal number of points to be transferred 
00144 !     directly to the destination process
00145 !----------------------------------------------------------------------------
00146 !
00147       ndir = 0
00148 !
00149          do ipart = 1, search%npart
00150          if (maxval(search%shape(:,:,ipart)-inter(:,:,ipart)) == 0) then
00151 !
00152 !           search%shape == inter; use intrinsic count function
00153 !
00154             ndir = ndir + count (search%search_mask(ipart)%vector(:))
00155          else
00156 !
00157 !           search%shape /= inter
00158 !
00159             call psmile_get_true_mask_entries (                 &
00160                         search%search_mask(ipart)%vector,       &
00161                         search%shape (:, :, ipart),             &
00162                         inter (:, :, ipart), n, ierror)
00163             if (ierror > 0) return
00164 
00165             ndir = ndir + n
00166          endif
00167 !
00168          end do ! ipart
00169 !
00170 !----------------------------------------------------------------------------
00171 !     Generate send areas for data exchange with application process
00172 !     ??? Derzeit werden die Daten Zellweise abgespeichert.
00173 !         Clustern der Daten, wenn moeglich !!!
00174 !----------------------------------------------------------------------------
00175 !
00176       if (ndir > 0) then
00177 !
00178         call psmile_get_info_index (method_id, send_direct_index, index, ierror)
00179         if (ierror > 0) return
00180 !
00181         dir_index = index
00182 !
00183         mp%send_infos_direct(index)%dest   = search%sender
00184         mp%send_infos_direct(index)%nvec   = 1
00185         mp%send_infos_direct(index)%nparts = 1
00186 !
00187         mp%send_infos_direct(index)%nrecv    = 0
00188         mp%send_infos_direct(index)%num2recv = 0
00189 !
00190         mp%send_infos_direct(index)%remote_method_id = search%msg_intersections%first_tgt_method_id
00191 !
00192 !       ... generate send areas for direct data exchange
00193 !       Hier nehme ich die search%search_mask statt Mask(mask_id)%mask_array
00194 !       weil die search mask kompakter (und damit schneller) sein sollte.
00195 !
00196         call psmile_locations_alloc (mp%send_infos_direct(index), ierror)
00197         if (ierror > 0) return
00198 !
00199         nprev = 0
00200 
00201            do ipart = 1, search%npart
00202            call psmile_store_dest_locs_3d_msk (                      &
00203                              search%search_mask(ipart)%vector,       &
00204                              search%shape (:, :, ipart),             &
00205                              inter (:, :, ipart),                    &
00206                              mp%send_infos_direct(index), ndir,      &
00207                              nprev, nadd, ierror)
00208            if (ierror > 0) return
00209 
00210            call psmile_store_source_locs_3d_msk (                    &
00211                              search%search_mask(ipart)%vector,       &
00212                              search%shape (:, :, ipart),             &
00213                              inter (:, :, ipart),                    &
00214                              mp%send_infos_direct(index), ndir,      &
00215                              nprev, nadd, ierror)
00216            if (ierror > 0) return
00217 !
00218            nprev = nprev + nadd
00219            end do ! ipart
00220 !
00221          mp%send_infos_direct(index)%nloc = nprev
00222 !
00223 !===> ... Transform indices in indices of target grid
00224 !
00225             do n = 1, ndim_3d
00226             mp%send_infos_direct(index)%dstijk (n, 1:nprev) = &
00227             mp%send_infos_direct(index)%dstijk (n, 1:nprev) - shift (n)
00228             end do
00229 !
00230       endif
00231 !
00232 !===> All done
00233 !
00234 #ifdef VERBOSE
00235       print 9980, trim(ch_id), ierror, dir_index
00236 
00237       call psmile_flushstd
00238 #endif /* VERBOSE */
00239 !
00240 !  Formats:
00241 !
00242 
00243 #ifdef VERBOSE
00244 
00245 9990 format (1x, a, ': psmile_locations_3d_mask: method_id', i3, &
00246                     ' to ', i3, '(', i2, ')')
00247 9980 format (1x, a, ': psmile_locations_3d_mask: eof ierror =', i3, &
00248                     '; dir_index', i7)
00249 
00250 #endif /* VERBOSE */
00251 
00252       end subroutine PSMILe_Locations_3d_mask

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1