psmile_search_donor_3d_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_Search_donor_3d_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_search_donor_3d_real (comp_info, &
00012                         found, locations, len, search,   &
00013                         field_list, n_vars,              &
00014                         grid_id, method_id, var_id, tol, ierror)
00015 !
00016 ! !USES:
00017 !
00018       use PRISM_constants
00019       use psmile_grid, only : common_grid_range
00020       use PSMILe, dummy_interface => PSMILe_Search_donor_3d_real
00021 
00022       implicit none
00023 !
00024 ! !INPUT PARAMETERS:
00025 !
00026       Type (Enddef_comp), Intent (In) :: comp_info
00027 
00028 !     Info on the component in which the donor cells
00029 !     should be searched.
00030 
00031       Type (Enddef_search)            :: search
00032 
00033 !     Info's on coordinates to be searched
00034 
00035       Integer, Intent (In)            :: len (search%npart)
00036 
00037 !     Number of coords to be searched
00038 !
00039       Integer, Intent (In)            :: grid_id
00040 !
00041 !     Grid id
00042 
00043       Integer, Intent (InOut)         :: method_id
00044 !
00045 !     Method id to be searched for.
00046 
00047       Integer, Intent (InOut)         :: var_id
00048 !
00049 !     Field Id (in this process) for the method
00050 !
00051       Integer, Intent (In)            :: n_vars
00052 !
00053 !     Number of additional fields to be searched for sending process
00054 !
00055       Integer, Intent (In)            :: field_list (nd_field_list, n_vars)
00056 !
00057 !     Info's on additional fields to be searched
00058 !   
00059       Real, Intent (In)               :: tol
00060 !
00061 !     Absolute tolerance for search of "identical" points
00062 !     TOL >= 0.0
00063 !
00064 ! !INPUT/OUTPUT PARAMETERS:
00065 !
00066       Type (integer_vector)           :: found (search%npart)
00067 
00068 !     Finest level number on which a grid cell was found for point I.
00069 !     Level number = -(nlev+1): Never found (input value)
00070 
00071       Type (integer_vector)           :: locations (search%npart)
00072 !
00073 !     Indices of the grid cell in which the point was found.
00074 !     Assumed input value locations(:)%vector (:, len) = 0
00075 !
00076 ! !OUTPUT PARAMETERS:
00077 !
00078       Integer, Intent (Out)           :: ierror
00079 
00080 !     Returns the error code of PSMILe_Search_donor_3d_real;
00081 !             ierror = 0 : No error
00082 !             ierror > 0 : Severe error
00083 !
00084 ! !LOCAL VARIABLES
00085 !
00086       Type (Corner_Block), Pointer :: corner_pointer
00087       Integer                      :: i, ipart
00088 !
00089 !     ... for locations searched
00090 !
00091       Integer                      :: ibeg, iend
00092       Integer                      :: control (2, ndim_3d, search%npart)
00093 !
00094       Type (real_vector)           :: array (ndim_3d, search%npart)
00095       Integer                      :: shape (2, ndim_3d, search%npart)
00096       Integer                      :: j, k, len1, len2, len3
00097 !
00098 !     ... for fields
00099 !
00100       Integer                      :: n_vars_ret
00101 !
00102 !     ... for methods
00103 !
00104       Integer                      :: dim
00105       Integer                      :: method_type, old
00106       Integer                      :: cpl_index, dir_index
00107       Integer                      :: len_cpl (search%npart)
00108       Integer                      :: m_levdim (ndim_3d)
00109       Type (Coords_Block),Pointer  :: coords_pointer
00110 !
00111       Logical                      :: changed
00112       Logical                      :: msk_required
00113       Type (integer_vector)        :: locations_save (search%npart)
00114       Type (integer_vector)        :: virtual_cell (search%npart) ! dummy
00115 !
00116       Type (Enddef_mg_real)      :: m_arrays
00117 !
00118       Real, Save                   :: period (ndim_3d)
00119 !
00120 !     ... for levels
00121 !
00122       Integer                      :: lev, levc, nlev, not_found
00123       Integer                      :: ijkinc (ndim_3d), ijkcoa (ndim_3d)
00124       Type(Grid), Pointer          :: grid_info
00125 !
00126 !     ... for error parameters
00127 !
00128       Integer, parameter           :: nerrp = 2
00129       Integer                      :: ierrp (nerrp)
00130 !
00131 ! !DESCRIPTION:
00132 !
00133 ! Subroutine "PSMILe_Search_donor_3d_real" searches the donor
00134 ! for the subgrid sent by the sending process.
00135 !
00136 ! Subroutine "PSMILe_Search_donor_irreg2_real" is designed for grids of
00137 ! type "PRISM_Irrlonlatvrt". We assume that this routine is only called
00138 ! for methods that are used for coupling.
00139 !
00140 !
00141 ! !REVISION HISTORY:
00142 !
00143 !   Date      Programmer   Description
00144 ! ----------  ----------   -----------
00145 ! 03.07.21    H. Ritzdorf  created
00146 !
00147 !EOP
00148 !----------------------------------------------------------------------
00149 !
00150 ! $Id: psmile_search_donor_3d_real.F90 2792 2010-12-01 17:35:05Z hanke $
00151 ! $Author: hanke $
00152 !
00153    Character(len=len_cvs_string), save :: mycvs = 
00154        '$Id: psmile_search_donor_3d_real.F90 2792 2010-12-01 17:35:05Z hanke $'
00155 !
00156 !----------------------------------------------------------------------
00157 !
00158 !  Initialization
00159 !
00160 #ifdef VERBOSE
00161       print 9990, trim(ch_id), comp_info%comp_id, search%sender
00162 
00163       call psmile_flushstd
00164 #endif /* VERBOSE */
00165 
00166       ierror = 0
00167       n_vars_ret = 0
00168       grid_info => Grids(grid_id)
00169       msk_required = .false.
00170 
00171       period (1:ndim_3d) = common_grid_range(2,1:ndim_3d) - &
00172                            common_grid_range(1,1:ndim_3d)
00173 !
00174       corner_pointer => Grids(grid_id)%corner_pointer
00175 !
00176 #ifdef PRISM_ASSERTION
00177       if (method_id > Number_of_Methods_allocated .or. &
00178           method_id < 1) then
00179          print *, trim(ch_id), "PSMILe_Search_donor_3d_real: method_id =", &
00180                   method_id
00181          call psmile_assert (__FILE__, __LINE__, &
00182               "Invalid Method id")
00183       endif
00184 !
00185       if ( grid_info%nlev <= 0 ) then
00186          call psmile_assert (__FILE__, __LINE__, "nlev <= 0")
00187       endif
00188 #endif
00189 !
00190 !-----------------------------------------------------------------------
00191 !     Generate 3d array if necessary
00192 !-----------------------------------------------------------------------
00193 !
00194       if (search%grid_type == PRISM_Reglonlatvrt .or. &
00195           search%grid_type == PRISM_Irrlonlat_regvrt) then
00196 !        ... Dimensions are only 1-dimensional or 2-dimensional
00197 !        ... Generate 3d arrays
00198 !
00199          shape = search%range(:, :, 1:search%npart)
00200 !
00201          do ipart = 1, search%npart
00202          Allocate (array(1,ipart)%vector(len(ipart)), STAT = ierror)
00203          if ( ierror > 0 ) then
00204             ierrp (1) = ierror
00205             ierrp (2) = len (ipart)
00206             ierror = PRISM_Error_Alloc
00207             call psmile_error ( ierror, 'array (1, ipart)', &
00208                                 ierrp, 2, __FILE__, __LINE__ )
00209             return
00210          endif
00211 !
00212          Allocate (array(2,ipart)%vector(len(ipart)), STAT = ierror)
00213             if ( ierror > 0 ) then
00214                ierrp (1) = ierror
00215             ierrp (2) = len (ipart)
00216                ierror = PRISM_Error_Alloc
00217             call psmile_error ( ierror, 'array2', &
00218                                    ierrp, 2, __FILE__, __LINE__ )
00219                return
00220             endif
00221 !
00222          Allocate (array(3,ipart)%vector(len(ipart)), STAT = ierror)
00223             if ( ierror > 0 ) then
00224                ierrp (1) = ierror
00225             ierrp (2) = len (ipart)
00226                ierror = PRISM_Error_Alloc
00227             call psmile_error ( ierror, 'array3', &
00228                                    ierrp, 2, __FILE__, __LINE__ )
00229                return
00230             endif
00231 !
00232 !===> ... Generate 3-dimensional temporary arrays
00233 !
00234          len1 =  search%range(2,1, ipart)-search%range(1,1, ipart) + 1
00235          len2 = (search%range(2,2, ipart)-search%range(1,2, ipart) + 1) * len1
00236 
00237          do k = 1, search%range(2,3,ipart)-search%range(1,3,ipart) + 1
00238             array(3,ipart)%vector((k-1)*len2+1:k*len2) = &
00239                search%search_real(3, ipart)%vector(k)
00240          end do
00241 
00242          if (search%grid_type == PRISM_Reglonlatvrt) then
00243 #ifdef PRISM_ASSERTION
00244             if (search%dims(1, ipart) /= len1) then
00245                call psmile_assert (__FILE__, __LINE__, &
00246                                    "dims(1, ipart) != len1")
00247             endif
00248 
00249             if (search%dims(1, ipart)*search%dims(2, ipart) /= len2) then
00250                call psmile_assert (__FILE__, __LINE__, &
00251                                    "dims(1, ipart)*dims(2, ipart) != len2")
00252             endif
00253 #endif
00254 !
00255             do k = 1, search%range(2,3, ipart)-search%range(1,3, ipart) + 1
00256                do j = 1, search%range(2,2, ipart)-search%range(1,2, ipart) + 1
00257                   array(1,ipart)%vector ((k-1)*len2+(j-1)*len1+1: &
00258                                          (k-1)*len2+j*len1) = &
00259                   search%search_real(1,ipart)%vector(1:len1)
00260                   array(2,ipart)%vector ((k-1)*len2+(j-1)*len1+1: &
00261                                          (k-1)*len2+j*len1) = &
00262                   search%search_real(2,ipart)%vector(j)
00263             end do
00264             end do
00265 
00266          else if (search%grid_type == PRISM_Irrlonlat_regvrt) then
00267 
00268 #ifdef PRISM_ASSERTION
00269             if (search%dims(1, ipart) /= len2) then
00270                print *, 'dims(1, ipart)', search%dims(1, ipart), len2
00271                call psmile_assert (__FILE__, __LINE__, &
00272                                    "dims(1, ipart) != len2")
00273             endif
00274 
00275             if (search%dims(2, ipart) /= len2) then
00276                print *, 'dims(2, ipart)', search%dims(2, ipart), len2
00277                call psmile_assert (__FILE__, __LINE__, &
00278                                    "dims(2, ipart) != len2")
00279             endif
00280 #endif
00281 
00282             do k = 1, search%range(2,3, ipart)-search%range(1,3, ipart) + 1
00283                do j = 1, len2
00284                   array(1,ipart)%vector ((k-1)*len2+1:k*len2) = &
00285                      search%search_real(1, ipart)%vector(1:len2)
00286                   array(2,ipart)%vector ((k-1)*len2+1:k*len2) = &
00287                      search%search_real(2, ipart)%vector(1:len2)
00288                end do
00289             end do
00290          endif
00291 !
00292          end do ! ipart
00293 !
00294       else
00295 !
00296             do ipart = 1, search%npart
00297             array(1,ipart)%vector => search%search_real(1,ipart)%vector
00298             array(2,ipart)%vector => search%search_real(2,ipart)%vector
00299             array(3,ipart)%vector => search%search_real(3,ipart)%vector
00300             end do ! ipart
00301 !
00302          shape = search%shape(1:2, 1:ndim_3d, 1:search%npart)
00303       endif
00304 !
00305 !-----------------------------------------------------------------------
00306 !     Coarsest level nlev
00307 !-----------------------------------------------------------------------
00308 !
00309       control (:, :, 1:search%npart) = search%range (:, :, 1:search%npart)
00310 !
00311       nlev = grid_info%nlev
00312 !
00313          do ipart = 1, search%npart
00314          lev = nlev
00315 !
00316          ibeg = 1
00317          iend = len (ipart)
00318 
00319 #ifdef PRISM_ASSERTION
00320          if (grid_info%mg_infos(lev)%levdim(1) /= 0 .or. &
00321              grid_info%mg_infos(lev)%levdim(2) /= 0 .or. &
00322              grid_info%mg_infos(lev)%levdim(3) /= 0) then
00323 
00324             call psmile_assert (__FILE__, __LINE__, &
00325                                 "coarsest level dim != 0")
00326          endif
00327 #endif
00328 !
00329          call psmile_mg_coarse_3d_real (lev, &
00330                      grid_info%mg_infos(lev)%real_arrays%chmin, &
00331                      grid_info%mg_infos(lev)%real_arrays%chmax, &
00332                      found(ipart)%vector, locations(ipart)%vector, &
00333                      array(1,ipart)%vector, array(2,ipart)%vector, &
00334                      array(3,ipart)%vector,  &
00335                      ibeg, iend)
00336 !
00337          if (ibeg > iend) then
00338 !           Make control range empty
00339             control (1, 1, ipart) = control (2, 1, ipart)
00340             control (2, 1, ipart) = control (2, 1, ipart) - 1
00341             cycle
00342          endif
00343 !
00344 !===> ... Find range
00345 !  den range besser einschraenken ! 3d Range !
00346 !
00347 !        call psmile_find_range (found(ipart)%vector, ibeg, iend,
00348 !                     control,
00349 !                     iret,  ierror)
00350 !        if (ierror > 0) return
00351 !
00352 !        if (iret == 0) go to 2000
00353 
00354 !
00355 !===> Multiple grids
00356 !
00357 !        ijkinc (1) = max (1, (iend-ibeg)/2)
00358 !        ijkinc (2) = max (1, (jend-jbeg)/2)
00359 !        ijkinc (3) = max (1, (kend-kbeg)/2)
00360 !
00361          ijkinc (:) = 1
00362 !
00363             do lev = nlev-1, 1, -1
00364 !
00365             levc = lev + 1
00366 !
00367             ijkcoa (:) = 2
00368 !
00369                do i = 1, ndim_3d
00370                if (grid_info%mg_infos(lev)%levdim(i) == &
00371                    grid_info%mg_infos(levc)%levdim(i)) then
00372                   ijkcoa (i) = 1
00373          endif
00374                enddo
00375 !
00376 ! chmin, chmax, midp als 4d array alloziieren !
00377 !
00378             call psmile_mg_next_level_3d_real ( grid_id, lev, nlev  ,  &
00379                 grid_info%mg_infos(lev)%real_arrays%chmin(1)%vector, &
00380                 grid_info%mg_infos(lev)%real_arrays%chmin(2)%vector, &
00381                 grid_info%mg_infos(lev)%real_arrays%chmin(3)%vector, &
00382                 grid_info%mg_infos(lev)%real_arrays%chmax(1)%vector, &
00383                 grid_info%mg_infos(lev)%real_arrays%chmax(2)%vector, &
00384                 grid_info%mg_infos(lev)%real_arrays%chmax(3)%vector, &
00385                 grid_info%mg_infos(lev)%real_arrays%midp(1)%vector,  &
00386                 grid_info%mg_infos(lev)%real_arrays%midp(2)%vector,  &
00387                 grid_info%mg_infos(lev)%real_arrays%midp(3)%vector,  &
00388                 grid_info%mg_infos(lev)%levdim,                        &
00389                 found(ipart)%vector, locations(ipart)%vector,          &
00390                 search%range(:,:, ipart),                              &
00391                 array(1,ipart)%vector, array(2,ipart)%vector,          &
00392                 array(3,ipart)%vector, shape (:, :, ipart),            &
00393                 control (:, :, ipart), ijkinc, ijkcoa, ierror)
00394 
00395             if (ierror > 0) return
00396 !
00397 !           ijkinc (:) = max (1, ijkinc(:)/ijkcoa(:))
00398             enddo ! lev
00399 !
00400 #ifdef SEARCH_ON_FINAL_GRID
00401 !
00402 !   Get final locations controlling the real cells
00403 !   Glaube nicht, dass das sich lohnt. Zu teuer. Lasse es trotzdem mal drin
00404 !
00405          if (corner_pointer%nbr_corners == 8) then
00406             call psmile_mg_final_level_3d_real ( comp_id, nlev,            &
00407                       found(ipart)%vector,  locations(ipart)%vector,       &
00408                       search%range (:, :, ipart),                          &
00409                       array(1,ipart)%vector, array(2,ipart)%vector,        &
00410                       array(3,ipart)%vector, shape (:, :, ipart), control, &
00411                       grid_id,                                             &
00412                       corner_pointer%corners_real(1)%vector,          &
00413                       corner_pointer%corners_real(2)%vector,          &
00414                       corner_pointer%corners_real(3)%vector,          &
00415                       corner_pointer%corner_shape,             &
00416                       corner_pointer%nbr_corners,              &
00417                       grid_info%ijk0, tol, ierror)
00418             if (ierror > 0) return
00419       endif
00420 #else
00421 !
00422 !   Compute locations
00423 !
00424 !cdir vector
00425           do i = 1, len (ipart)
00426           locations(ipart)%vector((i-1)*ndim_3d+1) = &
00427           locations(ipart)%vector((i-1)*ndim_3d+1) + grid_info%ijk0 (1)
00428           locations(ipart)%vector((i-1)*ndim_3d+2) = &
00429           locations(ipart)%vector((i-1)*ndim_3d+2) + grid_info%ijk0 (2)
00430           locations(ipart)%vector((i-1)*ndim_3d+3) = &
00431           locations(ipart)%vector((i-1)*ndim_3d+3) + grid_info%ijk0 (3)
00432          end do
00433 !
00434 #endif
00435          end do ! ipart
00436 !
00437 !-----------------------------------------------------------------------
00438 !     Control locations of method/field values
00439 !-----------------------------------------------------------------------
00440 !            
00441 !     changed = .true. : "found" and "locations" were changed by routine
00442 !                        PSMILe_mg_method_3d_real
00443 !
00444 !
00445       changed = .false.
00446 !
00447 !    ... Save locations since they may be changed by psmile_mg_method
00448 !        ??? Remove Points which are not found ???
00449 !        TODO: Abfrage auf n_vars verbessern
00450 !
00451       if (n_vars > 1) then
00452             do ipart = 1, search%npart
00453             Allocate (locations_save(ipart)%vector(1:len(ipart)), STAT = ierror)
00454             if ( ierror > 0 ) then
00455                ierrp (1) = ierror
00456             ierrp (2) = len(ipart)
00457                ierror = PRISM_Error_Alloc
00458             call psmile_error ( ierror, 'locations_save(ipart)%vector', &
00459                                    ierrp, 2, __FILE__, __LINE__ )
00460                return
00461             endif
00462             enddo 
00463 !
00464          do ipart = 1, search%npart
00465          locations_save(ipart)%vector(1:len(ipart)) = &
00466               locations(ipart)%vector(1:len(ipart))
00467          enddo
00468       endif
00469 !
00470 !     ... Get method type
00471 !
00472 1000  continue
00473 !
00474       method_type = Methods(method_id)%method_type
00475       coords_pointer => Methods(method_id)%coords_pointer
00476 
00477       if (method_type == PSMILe_PointMethod) then
00478 
00479          if (coords_pointer%coords_datatype /= MPI_REAL) then
00480             ierror = PRISM_Error_Internal
00481 
00482             call psmile_error ( ierror, 'Different datatypes in Grid and Method are currently not supported', &
00483                  ierrp, 0, __FILE__, __LINE__ )
00484          endif
00485 
00486 #ifdef PRISM_ASSERTION
00487             if (.not. Associated(coords_pointer%coords_real(1)%vector) ) then
00488             call psmile_assert (__FILE__, __LINE__, &
00489                  "x coordinates of method are not defined")
00490          endif
00491 #endif
00492 !
00493 !     ... reset found array
00494 !         ? changed parameter von psmile_mg_method, um das nur zu machen,
00495 !           wenn wert geaendert wurde ?
00496 !
00497          if (changed) then
00498                do ipart = 1, search%npart
00499                   found(ipart)%vector(1:len(ipart)) = abs (found(ipart)%vector(1:len(ipart)))
00500                enddo 
00501 !
00502                do ipart = 1, search%npart
00503                   locations     (ipart)%vector(1:len(ipart)) = &
00504                   locations_save(ipart)%vector(1:len(ipart))
00505                enddo
00506             endif
00507 !
00508 !     ... Compute min/max values of method grid
00509 !         Note: chmin, chmax and midp are dimensioned as
00510 !               grid_shape(1,*)-1:grid_shape(2,*)
00511 !         ??? Wirklich fuer das ganze Gitter berechnen oder
00512 !             Nur fuer das entsprechende Teilgitter unter Verwendung von
00513 !             grid_info%mg_infos(lev)%real_arrays%chmin(i)%vector,
00514 !             grid_info%mg_infos(lev)%real_arrays%chmax(i)%vector,
00515 !             grid_info%mg_infos(lev)%levdim(i)
00516 !    
00517          m_levdim (1:ndim_3d) = Grids(grid_id)%grid_shape(2,1:ndim_3d) - &
00518                                 Grids(grid_id)%grid_shape(1,1:ndim_3d) + 2
00519 
00520 ! 
00521          dim = m_levdim (1) * m_levdim (2) * m_levdim(3)
00522 
00523 !  ... Allocate arrays
00524 !
00525                do i = 1, ndim_3d
00526 
00527                Allocate (m_arrays%chmin(i)%vector(dim), STAT = ierror)
00528                if ( ierror > 0 ) then
00529                   ierrp (1) = ierror
00530                   ierrp (2) = dim
00531                   ierror = PRISM_Error_Alloc
00532                   call psmile_error ( ierror, 'm_arrays%chmin(i)%vector', &
00533                                       ierrp, 2, __FILE__, __LINE__ )
00534                   return
00535                endif
00536 
00537                Allocate (m_arrays%chmax(i)%vector(dim), STAT = ierror)
00538                if ( ierror > 0 ) then
00539                   ierrp (1) = ierror
00540                   ierrp (2) = dim
00541                   ierror = PRISM_Error_Alloc
00542                   call psmile_error ( ierror, 'm_arrays%chmax(i)%vector', &
00543                                       ierrp, 2, __FILE__, __LINE__ )
00544                   return
00545                endif
00546 
00547                Allocate (m_arrays%midp(i)%vector(dim), STAT = ierror)
00548                if ( ierror > 0 ) then
00549                   ierrp (1) = ierror
00550                   ierrp (2) = dim
00551                   ierror = PRISM_Error_Alloc
00552                   call psmile_error ( ierror, 'm_arrays%midp(i)%vector', &
00553                                       ierrp, 2, __FILE__, __LINE__ )
00554                   return
00555                endif
00556                end do
00557 !
00558 !           ... Compute bounding boxes for the cells
00559 !
00560                do i = 1, ndim_3d
00561                call psmile_bbcells_3d_real ( method_id,                    &
00562                   coords_pointer%coords_real(i)%vector,                    &
00563                   coords_pointer%coords_shape,                             &
00564                   Grids(grid_id)%grid_shape,                               &
00565                   corner_pointer%corners_real(i)%vector,                   &
00566                   corner_pointer%corner_shape(:,1:ndim_3d),                &
00567                   corner_pointer%nbr_corners,                              &
00568                   m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector,      &
00569                   m_arrays%midp(i)%vector,                                 &
00570                   m_levdim, grid_info%cyclic(i), period(i), ierror)
00571                end do
00572 !
00573 !           ... Search for coordinates in method grid
00574 !
00575                do ipart = 1, search%npart
00576                call psmile_mg_method_3d_real ( comp_info, nlev,       &
00577                   found(ipart)%vector,  locations(ipart)%vector,      &
00578                   search%range (:, :, ipart),                         &
00579                   array(1,ipart)%vector, array(2,ipart)%vector,       &
00580                   array(3,ipart)%vector, shape(:, :, ipart),          &
00581                   control (:, :, ipart),                              &
00582                   method_id,                                          &
00583                   coords_pointer%coords_real(1)%vector,               &
00584                   coords_pointer%coords_real(2)%vector,               &
00585                   coords_pointer%coords_real(3)%vector,               &
00586                   coords_pointer%coords_shape,                        &
00587                   grid_info%grid_shape,                               &
00588                   grid_info%cyclic,                                   &
00589                   m_arrays%chmin(1)%vector, m_arrays%chmin(2)%vector, &
00590                   m_arrays%chmin(3)%vector,                           &
00591                   m_arrays%chmax(1)%vector, m_arrays%chmax(2)%vector, &
00592                   m_arrays%chmax(3)%vector,                           &
00593                   m_arrays%midp (1)%vector, m_arrays%midp (2)%vector, &
00594                   m_arrays%midp (3)%vector,                           &
00595                   tol, ierror)
00596                if (ierror > 0) return
00597                end do ! ipart
00598 !
00599             changed = .true.
00600 !
00601 !           ... Free arrays allocated
00602 !
00603                do i = 1, ndim_3d
00604                Deallocate (m_arrays%midp (i)%vector)
00605                Deallocate (m_arrays%chmin(i)%vector)
00606                Deallocate (m_arrays%chmax(i)%vector)
00607                end do
00608 
00609       else if (method_type == PSMILe_VectorPointMethod) then
00610 
00611          ierrp (1) = method_type
00612          ierror = PRISM_Error_Internal
00613          call psmile_error ( ierror, 'Vector Method is currently not supported', &
00614               ierrp, 1, __FILE__, __LINE__ )
00615 
00616       else if (method_type == PSMILe_SubgridMethod) then
00617          ierrp (1) = method_type
00618          ierror = PRISM_Error_Internal
00619          call psmile_error ( ierror, 'Subgrid Method is currently not supported', &
00620               ierrp, 1, __FILE__, __LINE__ )
00621 
00622       else
00623          ierrp (1) = method_type
00624          ierror = PRISM_Error_Internal
00625          call psmile_error ( ierror, 'unsupported method type', &
00626               ierrp, 1, __FILE__, __LINE__ )
00627 
00628       endif
00629 !
00630 !===> Send information on locations found back to sending process
00631 !
00632 ! Hier gibt es verschieden Moeglichkeiten. Derzeit ist mir nicht klar,
00633 ! was man machen muss (das haengt auch davon ab, was in PRISM erlaubt ist).
00634 !
00635 ! (i)  Die Applikationen einigen sich zuerst auf die Punkte,
00636 !      welche Werte von wem geupdated werden, wenn es mehrere Moeglich-
00637 !      keiten gibt, und der empfangende Process teilt dem sendenden
00638 !      Process mit, welche Werte tatsaechlich gesendet werden sollen.
00639 !
00640 ! (ii) Man sendet die Informationen an den Coupler und der Coupler
00641 !      ist dafuer zustaendig.
00642 !      Das klaert noch nicht das Problem, dass ohne
00643 !      Verwendung des Coupler Werte von unterschiedlichen Applikationen
00644 !      geupdated werden koennen.
00645 !
00646 ! Da muss man sich mal drueber unterhalten oder dass muss sich mit den
00647 ! Applikationen herausstellen. Zunaechts wird (ii) realisiert.
00648 !
00649 !
00650 1500  continue
00651 !
00652 !     Eliminate doubly defined entries
00653 !
00654          call psmile_compact_locations ( grid_id, search, ndim_1d, found, ierror )
00655 !
00656 !     Apply mask on points to be searched if required;
00657 !     i.e. set "found" to -(nlev+1) ("not found")
00658 !          if search mask is not true.
00659 !
00660          if (Associated(search%search_mask)) then
00661 
00662 #ifdef PRISM_ASSERTION
00663 !     Vorsicht: Hier wird search%shape == search_range angenommen.
00664             if (maxval (search%shape-search%range) /= 0) then
00665                call psmile_assert (__FILE__, __LINE__, &
00666                                    "search%shape != search%range")
00667             endif
00668 #endif
00669             not_found = - (grid_info%nlev+1)
00670 
00671             do ipart = 1, search%npart
00672             len3 = (search%range(2,1,ipart)-search%range(1,1,ipart)+1) * &
00673                    (search%range(2,2,ipart)-search%range(1,2,ipart)+1) * &
00674                    (search%range(2,3,ipart)-search%range(1,3,ipart)+1)
00675 !cdir vector
00676                do i = 1, len3
00677                if (.not. search%search_mask(ipart)%vector(i)) then
00678                   found(ipart)%vector(i) = not_found
00679                endif
00680                end do ! i
00681             end do ! ipart
00682          endif
00683 !
00684 !        Store source and destination locations on data found
00685 !
00686          call psmile_locations_3d (found, locations, search%range, &
00687                    control, search, method_id, msk_required,       &
00688                    virtual_cell, .false.,                          &
00689                    dir_index, cpl_index, len_cpl, ierror)
00690          if (ierror > 0) return
00691 
00692          if (cpl_index > 0) then
00693 
00694 !              Transfer info to the coupler process
00695             
00696                call psmile_info_trs_locs_3d_real (comp_info, &
00697                       array, shape, control, len_cpl,        &
00698                  var_id, Grids(grid_id)%grid_shape,          &
00699                  search, method_id,                          &
00700                  cpl_index, ierror)
00701             if (ierror > 0) return
00702 !
00703          endif
00704 !
00705 !        Store send info's for field "var_id"
00706 !
00707          call psmile_store_send_info (var_id, search%msg_intersections%transient_out_id, &
00708                                       dir_index, cpl_index, PRISM_UNDEFINED, ierror)
00709          if (ierror > 0) return
00710 !
00711 !        Send locations to the destination process
00712 !
00713          call psmile_return_locations_3d (search%msg_intersections,     &
00714                                           search%sender, method_id,     &
00715                                           dir_index, cpl_index, n_vars, &
00716                                           n_vars_ret, ierror)
00717          if (ierror > 0) return
00718 !
00719 !===>    Get next field/var id and mask id for the target process.
00720 !        Currently we DON'T look for the best field
00721 !
00722          if (n_vars_ret < n_vars) then
00723 !
00724             call psmile_get_next_field (comp_info, search, field_list, &
00725                                         n_vars, n_vars_ret, var_id, ierror)
00726             if (ierror > 0) return
00727 !
00728             old = method_id
00729             method_id = Fields(var_id)%method_id
00730 !           grid_id = Methods(method_id)%grid_id
00731 !           datatype = Grids(grid_id)%corner_pointer%corner_datatype
00732 !
00733             if (old == method_id) go to 1500
00734             go to 1000
00735          endif
00736 !
00737 !    ... Free locations saved
00738 !
00739       if (n_vars > 1) then
00740          do ipart = 1, search%npart
00741          Deallocate (locations_save(ipart)%vector)
00742          enddo
00743       endif
00744 !
00745 !-----------------------------------------------------------------------
00746 !     All done; Free memory allocated
00747 !-----------------------------------------------------------------------
00748 !
00749       if (search%grid_type == PRISM_Reglonlatvrt .or. &
00750           search%grid_type == PRISM_Irrlonlat_regvrt) then
00751             do ipart = 1, search%npart
00752             Deallocate (array(3,ipart)%vector)
00753             Deallocate (array(2,ipart)%vector)
00754             Deallocate (array(1,ipart)%vector)
00755          end do ! ipart
00756       endif
00757 !
00758 #ifdef VERBOSE
00759       print 9980, trim(ch_id), grid_info%comp_id, search%sender, ierror
00760 
00761       call psmile_flushstd
00762 #endif /* VERBOSE */
00763 !
00764 !  Formats:
00765 !
00766 9990 format (1x, a, ': psmile_search_donor_3d_real: comp_id =', i3, &
00767                     '; sender =', i4)
00768 9980 format (1x, a, ': psmile_search_donor_3d_real: comp_id =', i3, &
00769                     '; eof sender =', i3, ', ierror =', i4)
00770 
00771       end subroutine PSMILe_Search_donor_3d_real

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1