psmile_search_donor_irreg3_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_irreg3_real
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_search_donor_irreg3_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_irreg3_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%search_data%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       Type (enddef_field_info), Intent (In) :: 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), Intent (InOut) :: found (search%search_data%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), Intent (InOut) :: locations (search%search_data%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_irreg3_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                      :: control (2, ndim_3d, search%search_data%npart)
00092 !
00093       Type (real_vector)           :: array (ndim_3d, search%search_data%npart)
00094       Integer                      :: shape (2, ndim_3d, search%search_data%npart)
00095       Integer                      :: j, k, len1, len2, len3
00096 !
00097 !     ... for fields
00098 !
00099       Integer                      :: n_vars_ret
00100 !
00101 !     ... for methods
00102 !
00103       Integer                      :: dim
00104       Integer                      :: method_type, old
00105       Integer                      :: cpl_index, dir_index
00106       Integer                      :: len_cpl (search%search_data%npart)
00107       Integer                      :: m_levdim (ndim_3d)
00108       Type (Coords_Block),Pointer  :: coords_pointer
00109 !
00110       Logical                      :: changed
00111       Logical                      :: msk_required
00112       Type (integer_vector)        :: locations_save (search%search_data%npart)
00113       Type (integer_vector)        :: virtual_cell (search%search_data%npart) ! dummy
00114 !
00115       Type (Enddef_mg_real)      :: m_arrays
00116 !
00117       Real, Save                   :: period (ndim_3d)
00118 !
00119 !     ... for levels
00120 !
00121       Integer                      :: nlev, not_found
00122       Type(Grid), Pointer          :: grid_info
00123 !
00124 !     ... for error parameters
00125 !
00126       Integer, parameter           :: nerrp = 2
00127       Integer                      :: ierrp (nerrp)
00128 !
00129 ! !DESCRIPTION:
00130 !
00131 ! Subroutine "psmile_search_donor_irreg3_real" searches the donor
00132 ! for the subgrid sent by the sending process.
00133 !
00134 ! Subroutine "PSMILe_Search_donor_irreg2_real" is designed for grids of
00135 ! type "PRISM_Irrlonlatvrt". We assume that this routine is only called
00136 ! for methods that are used for coupling.
00137 !
00138 !
00139 ! !REVISION HISTORY:
00140 !
00141 !   Date      Programmer   Description
00142 ! ----------  ----------   -----------
00143 ! 03.07.21    H. Ritzdorf  created as psmile_search_donor_3d_dble.F90
00144 ! 20.04.2011  M. Hanke     extracted from psmile_search_donor_3d_dble.F90
00145 !
00146 !EOP
00147 !----------------------------------------------------------------------
00148 !
00149 ! $Id: psmile_search_donor_irreg3_real.F90 3174 2011-05-04 09:41:43Z hanke $
00150 ! $Author: hanke $
00151 !
00152    Character(len=len_cvs_string), save :: mycvs = 
00153        '$Id: psmile_search_donor_irreg3_real.F90 3174 2011-05-04 09:41:43Z hanke $'
00154 !
00155 !----------------------------------------------------------------------
00156 !
00157 !  Initialization
00158 !
00159 #ifdef VERBOSE
00160       print 9990, trim(ch_id), comp_info%comp_id, search%sender
00161 
00162       call psmile_flushstd
00163 #endif /* VERBOSE */
00164 
00165       ierror = 0
00166       n_vars_ret = 0
00167       grid_info => Grids(grid_id)
00168       msk_required = .false.
00169 
00170       period (1:ndim_3d) = common_grid_range(2,1:ndim_3d) - &
00171                            common_grid_range(1,1:ndim_3d)
00172 !
00173       corner_pointer => Grids(grid_id)%corner_pointer
00174 !
00175 #ifdef PRISM_ASSERTION
00176       if (method_id > Number_of_Methods_allocated .or. &
00177           method_id < 1) then
00178          print *, trim(ch_id), "psmile_search_donor_irreg3_real: method_id =", &
00179                   method_id
00180          call psmile_assert (__FILE__, __LINE__, &
00181               "Invalid Method id")
00182       endif
00183 !
00184       if ( grid_info%nlev <= 0 ) then
00185          call psmile_assert (__FILE__, __LINE__, "nlev <= 0")
00186       endif
00187 #endif
00188 !
00189 !-----------------------------------------------------------------------
00190 !     Generate 3d array if necessary
00191 !-----------------------------------------------------------------------
00192 !
00193       if (search%search_data%grid_type == PRISM_Reglonlatvrt .or. &
00194           search%search_data%grid_type == PRISM_Irrlonlat_regvrt) then
00195 !        ... Dimensions are only 1-dimensional or 2-dimensional
00196 !        ... Generate 3d arrays
00197 !
00198          shape = search%search_data%range(:, :, 1:search%search_data%npart)
00199 !
00200          do ipart = 1, search%search_data%npart
00201             Allocate (array(1,ipart)%vector(len(ipart)), STAT = ierror)
00202             if ( ierror > 0 ) then
00203                ierrp (1) = ierror
00204                ierrp (2) = len (ipart)
00205                ierror = PRISM_Error_Alloc
00206                call psmile_error ( ierror, 'array (1, ipart)', &
00207                                    ierrp, 2, __FILE__, __LINE__ )
00208                return
00209             endif
00210 !
00211             Allocate (array(2,ipart)%vector(len(ipart)), STAT = ierror)
00212             if ( ierror > 0 ) then
00213                ierrp (1) = ierror
00214             ierrp (2) = len (ipart)
00215                ierror = PRISM_Error_Alloc
00216             call psmile_error ( ierror, 'array2', &
00217                                    ierrp, 2, __FILE__, __LINE__ )
00218                return
00219             endif
00220 !
00221             Allocate (array(3,ipart)%vector(len(ipart)), STAT = ierror)
00222             if ( ierror > 0 ) then
00223                ierrp (1) = ierror
00224             ierrp (2) = len (ipart)
00225                ierror = PRISM_Error_Alloc
00226             call psmile_error ( ierror, 'array3', &
00227                                    ierrp, 2, __FILE__, __LINE__ )
00228                return
00229             endif
00230 !
00231 !===> ... Generate 3-dimensional temporary arrays
00232 !
00233             len1 =  search%search_data%range(2,1, ipart)-search%search_data%range(1,1, ipart) + 1
00234             len2 = (search%search_data%range(2,2, ipart)-search%search_data%range(1,2, ipart) + 1) * len1
00235 
00236             do k = 1, search%search_data%range(2,3,ipart)-search%search_data%range(1,3,ipart) + 1
00237                array(3,ipart)%vector((k-1)*len2+1:k*len2) = &
00238                   search%search_data%search_real(3, ipart)%vector(k)
00239             end do
00240 
00241             if (search%search_data%grid_type == PRISM_Reglonlatvrt) then
00242 #ifdef PRISM_ASSERTION
00243                if (search%search_data%dim_size(1, ipart) /= len1) then
00244                   call psmile_assert (__FILE__, __LINE__, &
00245                                       "dims(1, ipart) != len1")
00246                endif
00247 
00248                if (search%search_data%dim_size(1, ipart)*search%search_data%dim_size(2, ipart) /= len2) then
00249                   call psmile_assert (__FILE__, __LINE__, &
00250                                       "dims(1, ipart)*dims(2, ipart) != len2")
00251                endif
00252 #endif
00253 !
00254                do k = 1, search%search_data%range(2,3, ipart)-search%search_data%range(1,3, ipart) + 1
00255                   do j = 1, search%search_data%range(2,2, ipart)-search%search_data%range(1,2, ipart) + 1
00256                      array(1,ipart)%vector ((k-1)*len2+(j-1)*len1+1: &
00257                                             (k-1)*len2+j*len1) = &
00258                      search%search_data%search_real(1,ipart)%vector(1:len1)
00259                      array(2,ipart)%vector ((k-1)*len2+(j-1)*len1+1: &
00260                                             (k-1)*len2+j*len1) = &
00261                      search%search_data%search_real(2,ipart)%vector(j)
00262                   end do
00263                end do
00264 
00265             else if (search%search_data%grid_type == PRISM_Irrlonlat_regvrt) then
00266 
00267 #ifdef PRISM_ASSERTION
00268                if (search%search_data%dim_size(1, ipart) /= len2) then
00269                   print *, 'dims(1, ipart)', search%search_data%dim_size(1, ipart), len2
00270                   call psmile_assert (__FILE__, __LINE__, &
00271                                       "dims(1, ipart) != len2")
00272                endif
00273 
00274                if (search%search_data%dim_size(2, ipart) /= len2) then
00275                   print *, 'dims(2, ipart)', search%search_data%dim_size(2, ipart), len2
00276                   call psmile_assert (__FILE__, __LINE__, &
00277                                       "dims(2, ipart) != len2")
00278                endif
00279 #endif
00280 
00281                do k = 1, search%search_data%range(2,3, ipart)-search%search_data%range(1,3, ipart) + 1
00282                   do j = 1, len2
00283                      array(1,ipart)%vector ((k-1)*len2+1:k*len2) = &
00284                         search%search_data%search_real(1, ipart)%vector(1:len2)
00285                      array(2,ipart)%vector ((k-1)*len2+1:k*len2) = &
00286                         search%search_data%search_real(2, ipart)%vector(1:len2)
00287                   end do
00288                end do
00289             endif
00290 !
00291          end do ! ipart
00292 !
00293       else
00294 !
00295          do ipart = 1, search%search_data%npart
00296             array(1,ipart)%vector => search%search_data%search_real(1,ipart)%vector
00297             array(2,ipart)%vector => search%search_data%search_real(2,ipart)%vector
00298             array(3,ipart)%vector => search%search_data%search_real(3,ipart)%vector
00299          end do ! ipart
00300 !
00301          shape = search%search_data%shape(1:2, 1:ndim_3d, 1:search%search_data%npart)
00302       endif
00303 !
00304 !-----------------------------------------------------------------------
00305 !     Control locations of method/field values
00306 !-----------------------------------------------------------------------
00307 !            
00308 !     changed = .true. : "found" and "locations" were changed by routine
00309 !                        PSMILe_mg_method_3d_real
00310 !
00311 !
00312       changed = .false.
00313 !
00314 !    ... Save locations since they may be changed by psmile_mg_method
00315 !        ??? Remove Points which are not found ???
00316 !        TODO: Abfrage auf n_vars verbessern
00317 !
00318       if (n_vars > 1) then
00319             do ipart = 1, search%search_data%npart
00320             Allocate (locations_save(ipart)%vector(1:len(ipart)), STAT = ierror)
00321             if ( ierror > 0 ) then
00322                ierrp (1) = ierror
00323             ierrp (2) = len(ipart)
00324                ierror = PRISM_Error_Alloc
00325             call psmile_error ( ierror, 'locations_save(ipart)%vector', &
00326                                    ierrp, 2, __FILE__, __LINE__ )
00327                return
00328             endif
00329             enddo 
00330 !
00331          do ipart = 1, search%search_data%npart
00332          locations_save(ipart)%vector(1:len(ipart)) = &
00333               locations(ipart)%vector(1:len(ipart))
00334          enddo
00335       endif
00336 !
00337 !     ... Get method type
00338 !
00339 1000  continue
00340 !
00341       method_type = Methods(method_id)%method_type
00342       coords_pointer => Methods(method_id)%coords_pointer
00343 
00344       if (method_type == PSMILe_PointMethod) then
00345 
00346          if (coords_pointer%coords_datatype /= MPI_REAL) then
00347             ierror = PRISM_Error_Internal
00348 
00349             call psmile_error ( ierror, 'Different datatypes in Grid and Method are currently not supported', &
00350                  ierrp, 0, __FILE__, __LINE__ )
00351          endif
00352 
00353 #ifdef PRISM_ASSERTION
00354             if (.not. Associated(coords_pointer%coords_real(1)%vector) ) then
00355             call psmile_assert (__FILE__, __LINE__, &
00356                  "x coordinates of method are not defined")
00357          endif
00358 #endif
00359 !
00360 !     ... reset found array
00361 !         ? changed parameter von psmile_mg_method, um das nur zu machen,
00362 !           wenn wert geaendert wurde ?
00363 !
00364          if (changed) then
00365                do ipart = 1, search%search_data%npart
00366                   found(ipart)%vector(1:len(ipart)) = abs (found(ipart)%vector(1:len(ipart)))
00367                enddo 
00368 !
00369                do ipart = 1, search%search_data%npart
00370                   locations     (ipart)%vector(1:len(ipart)) = &
00371                   locations_save(ipart)%vector(1:len(ipart))
00372                enddo
00373             endif
00374 !
00375 !     ... Compute min/max values of method grid
00376 !         Note: chmin, chmax and midp are dimensioned as
00377 !               grid_shape(1,*)-1:grid_shape(2,*)
00378 !         ??? Wirklich fuer das ganze Gitter berechnen oder
00379 !             Nur fuer das entsprechende Teilgitter unter Verwendung von
00380 !             grid_info%mg_infos(lev)%real_arrays%chmin(i)%vector,
00381 !             grid_info%mg_infos(lev)%real_arrays%chmax(i)%vector,
00382 !             grid_info%mg_infos(lev)%levdim(i)
00383 !    
00384          m_levdim (1:ndim_3d) = Grids(grid_id)%grid_shape(2,1:ndim_3d) - &
00385                                 Grids(grid_id)%grid_shape(1,1:ndim_3d) + 2
00386 
00387 ! 
00388          dim = m_levdim (1) * m_levdim (2) * m_levdim(3)
00389 
00390 !  ... Allocate arrays
00391 !
00392                do i = 1, ndim_3d
00393 
00394                Allocate (m_arrays%chmin(i)%vector(dim), STAT = ierror)
00395                if ( ierror > 0 ) then
00396                   ierrp (1) = ierror
00397                   ierrp (2) = dim
00398                   ierror = PRISM_Error_Alloc
00399                   call psmile_error ( ierror, 'm_arrays%chmin(i)%vector', &
00400                                       ierrp, 2, __FILE__, __LINE__ )
00401                   return
00402                endif
00403 
00404                Allocate (m_arrays%chmax(i)%vector(dim), STAT = ierror)
00405                if ( ierror > 0 ) then
00406                   ierrp (1) = ierror
00407                   ierrp (2) = dim
00408                   ierror = PRISM_Error_Alloc
00409                   call psmile_error ( ierror, 'm_arrays%chmax(i)%vector', &
00410                                       ierrp, 2, __FILE__, __LINE__ )
00411                   return
00412                endif
00413 
00414                Allocate (m_arrays%midp(i)%vector(dim), STAT = ierror)
00415                if ( ierror > 0 ) then
00416                   ierrp (1) = ierror
00417                   ierrp (2) = dim
00418                   ierror = PRISM_Error_Alloc
00419                   call psmile_error ( ierror, 'm_arrays%midp(i)%vector', &
00420                                       ierrp, 2, __FILE__, __LINE__ )
00421                   return
00422                endif
00423                end do
00424 !
00425 !           ... Compute bounding boxes for the cells
00426 !
00427                do i = 1, ndim_3d
00428                call psmile_bbcells_3d_real ( method_id,                    &
00429                   coords_pointer%coords_real(i)%vector,                    &
00430                   coords_pointer%coords_shape,                             &
00431                   Grids(grid_id)%grid_shape,                               &
00432                   corner_pointer%corners_real(i)%vector,                   &
00433                   corner_pointer%corner_shape(:,1:ndim_3d),                &
00434                   corner_pointer%nbr_corners,                              &
00435                   m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector,      &
00436                   m_arrays%midp(i)%vector,                                 &
00437                   m_levdim, grid_info%cyclic(i), period(i), ierror)
00438                end do
00439 !
00440 !           ... Search for coordinates in method grid
00441 !
00442                do ipart = 1, search%search_data%npart
00443                call psmile_mg_method_3d_real ( comp_info, nlev,       &
00444                   found(ipart)%vector,  locations(ipart)%vector,      &
00445                   search%search_data%range (:, :, ipart),                         &
00446                   array(1,ipart)%vector, array(2,ipart)%vector,       &
00447                   array(3,ipart)%vector, shape(:, :, ipart),          &
00448                   control (:, :, ipart),                              &
00449                   method_id,                                          &
00450                   coords_pointer%coords_real(1)%vector,               &
00451                   coords_pointer%coords_real(2)%vector,               &
00452                   coords_pointer%coords_real(3)%vector,               &
00453                   coords_pointer%coords_shape,                        &
00454                   grid_info%grid_shape,                               &
00455                   grid_info%cyclic,                                   &
00456                   m_arrays%chmin(1)%vector, m_arrays%chmin(2)%vector, &
00457                   m_arrays%chmin(3)%vector,                           &
00458                   m_arrays%chmax(1)%vector, m_arrays%chmax(2)%vector, &
00459                   m_arrays%chmax(3)%vector,                           &
00460                   m_arrays%midp (1)%vector, m_arrays%midp (2)%vector, &
00461                   m_arrays%midp (3)%vector,                           &
00462                   tol, ierror)
00463                if (ierror > 0) return
00464                end do ! ipart
00465 !
00466             changed = .true.
00467 !
00468 !           ... Free arrays allocated
00469 !
00470                do i = 1, ndim_3d
00471                Deallocate (m_arrays%midp (i)%vector)
00472                Deallocate (m_arrays%chmin(i)%vector)
00473                Deallocate (m_arrays%chmax(i)%vector)
00474                end do
00475 
00476       else if (method_type == PSMILe_VectorPointMethod) then
00477 
00478          ierrp (1) = method_type
00479          ierror = PRISM_Error_Internal
00480          call psmile_error ( ierror, 'Vector Method is currently not supported', &
00481               ierrp, 1, __FILE__, __LINE__ )
00482 
00483       else if (method_type == PSMILe_SubgridMethod) then
00484          ierrp (1) = method_type
00485          ierror = PRISM_Error_Internal
00486          call psmile_error ( ierror, 'Subgrid Method is currently not supported', &
00487               ierrp, 1, __FILE__, __LINE__ )
00488 
00489       else
00490          ierrp (1) = method_type
00491          ierror = PRISM_Error_Internal
00492          call psmile_error ( ierror, 'unsupported method type', &
00493               ierrp, 1, __FILE__, __LINE__ )
00494 
00495       endif
00496 !
00497 !===> Send information on locations found back to sending process
00498 !
00499 ! Hier gibt es verschieden Moeglichkeiten. Derzeit ist mir nicht klar,
00500 ! was man machen muss (das haengt auch davon ab, was in PRISM erlaubt ist).
00501 !
00502 ! (i)  Die Applikationen einigen sich zuerst auf die Punkte,
00503 !      welche Werte von wem geupdated werden, wenn es mehrere Moeglich-
00504 !      keiten gibt, und der empfangende Process teilt dem sendenden
00505 !      Process mit, welche Werte tatsaechlich gesendet werden sollen.
00506 !
00507 ! (ii) Man sendet die Informationen an den Coupler und der Coupler
00508 !      ist dafuer zustaendig.
00509 !      Das klaert noch nicht das Problem, dass ohne
00510 !      Verwendung des Coupler Werte von unterschiedlichen Applikationen
00511 !      geupdated werden koennen.
00512 !
00513 ! Da muss man sich mal drueber unterhalten oder dass muss sich mit den
00514 ! Applikationen herausstellen. Zunaechts wird (ii) realisiert.
00515 !
00516 !
00517 1500  continue
00518 !
00519 !     Eliminate doubly defined entries
00520 !
00521          call psmile_compact_locations ( grid_id, search, ndim_1d, found, ierror )
00522 !
00523 !     Apply mask on points to be searched if required;
00524 !     i.e. set "found" to -(nlev+1) ("not found")
00525 !          if search mask is not true.
00526 !
00527          if (Associated(search%search_mask)) then
00528 
00529 #ifdef PRISM_ASSERTION
00530 !     Vorsicht: Hier wird search%search_data%shape == search_range angenommen.
00531             if (maxval (search%search_data%shape-search%search_data%range) /= 0) then
00532                call psmile_assert (__FILE__, __LINE__, &
00533                                    "search%search_data%shape != search%search_data%range")
00534             endif
00535 #endif
00536             not_found = - (grid_info%nlev+1)
00537 
00538             do ipart = 1, search%search_data%npart
00539             len3 = (search%search_data%range(2,1,ipart)-search%search_data%range(1,1,ipart)+1) * &
00540                    (search%search_data%range(2,2,ipart)-search%search_data%range(1,2,ipart)+1) * &
00541                    (search%search_data%range(2,3,ipart)-search%search_data%range(1,3,ipart)+1)
00542 !cdir vector
00543                do i = 1, len3
00544                if (.not. search%search_mask(ipart)%vector(i)) then
00545                   found(ipart)%vector(i) = not_found
00546                endif
00547                end do ! i
00548             end do ! ipart
00549          endif
00550 !
00551 !        Store source and destination locations on data found
00552 !
00553          call psmile_locations_3d (found, locations, search%search_data%range, &
00554                    control, search, method_id, msk_required,       &
00555                    virtual_cell, .false.,                          &
00556                    dir_index, cpl_index, len_cpl, ierror)
00557          if (ierror > 0) return
00558 
00559          if (cpl_index > 0) then
00560 
00561 !              Transfer info to the coupler process
00562             
00563                call psmile_info_trs_locs_3d_real (comp_info, &
00564                       array, shape, control, len_cpl,        &
00565                  var_id, Grids(grid_id)%grid_shape,          &
00566                  search, method_id,                          &
00567                  cpl_index, ierror)
00568             if (ierror > 0) return
00569 !
00570          endif
00571 !
00572 !        Store send info's for field "var_id"
00573 !
00574          call psmile_store_send_info (var_id, search%msg_intersections%field_info%transient_out_id, &
00575                                       dir_index, cpl_index, PRISM_UNDEFINED, ierror)
00576          if (ierror > 0) return
00577 !
00578 !        Send locations to the destination process
00579 !
00580          call psmile_return_locations_3d (search%msg_intersections,     &
00581                                           search%sender, method_id,     &
00582                                           dir_index, cpl_index, n_vars, &
00583                                           n_vars_ret, ierror)
00584          if (ierror > 0) return
00585 !
00586 !===>    Get next field/var id and mask id for the target process.
00587 !        Currently we DON'T look for the best field
00588 !
00589          if (n_vars_ret < n_vars) then
00590 !
00591             call psmile_get_next_field (comp_info, search, field_list, &
00592                                         n_vars, n_vars_ret, var_id, ierror)
00593             if (ierror > 0) return
00594 !
00595             old = method_id
00596             method_id = Fields(var_id)%method_id
00597 !           grid_id = Methods(method_id)%grid_id
00598 !           datatype = Grids(grid_id)%corner_pointer%corner_datatype
00599 !
00600             if (old == method_id) go to 1500
00601             go to 1000
00602          endif
00603 !
00604 !    ... Free locations saved
00605 !
00606       if (n_vars > 1) then
00607          do ipart = 1, search%search_data%npart
00608          Deallocate (locations_save(ipart)%vector)
00609          enddo
00610       endif
00611 !
00612 !-----------------------------------------------------------------------
00613 !     All done; Free memory allocated
00614 !-----------------------------------------------------------------------
00615 !
00616       if (search%search_data%grid_type == PRISM_Reglonlatvrt .or. &
00617           search%search_data%grid_type == PRISM_Irrlonlat_regvrt) then
00618             do ipart = 1, search%search_data%npart
00619             Deallocate (array(3,ipart)%vector)
00620             Deallocate (array(2,ipart)%vector)
00621             Deallocate (array(1,ipart)%vector)
00622          end do ! ipart
00623       endif
00624 !
00625 #ifdef VERBOSE
00626       print 9980, trim(ch_id), grid_info%comp_id, search%sender, ierror
00627 
00628       call psmile_flushstd
00629 #endif /* VERBOSE */
00630 !
00631 !  Formats:
00632 !
00633 9990 format (1x, a, ': psmile_search_donor_irreg3_real: comp_id =', i3, &
00634                     '; sender =', i4)
00635 9980 format (1x, a, ': psmile_search_donor_irreg3_real: comp_id =', i3, &
00636                     '; eof sender =', i3, ', ierror =', i4)
00637 
00638       end subroutine psmile_search_donor_irreg3_real

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1