00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011       subroutine psmile_search_donor_gridless (comp_info, search, &
00012                               field_list, n_vars,                 &
00013                               grid_id, method_id, var_id, ierror)
00014 
00015 
00016 
00017       use PRISM_constants
00018 
00019       use PSMILe, dummy_interface => PSMILe_Search_donor_gridless
00020 
00021       implicit none
00022 
00023 
00024 
00025       Type (Enddef_comp),   Intent (In)    :: comp_info
00026 
00027 
00028 
00029       Integer,              Intent (In)    :: grid_id
00030 
00031 
00032 
00033       Integer,              Intent (In)    :: n_vars
00034 
00035 
00036 
00037       Integer,              Intent (In)    :: field_list (nd_field_list, n_vars)
00038 
00039 
00040 
00041 
00042 
00043       Type (Enddef_search), Intent (InOut) :: search
00044 
00045 
00046 
00047       Integer,              Intent (InOut) :: method_id
00048 
00049 
00050 
00051       Integer,              Intent (InOut) :: var_id
00052 
00053 
00054 
00055 
00056 
00057       Integer,              Intent (Out)   :: ierror
00058 
00059 
00060 
00061 
00062 
00063 
00064 
00065       Integer, Parameter              :: cpl_index = PRISM_UNDEFINED
00066 
00067 
00068 
00069       Integer                         :: n, nb, nl, nbr_blocks
00070       Integer                         :: i, ipart, npart, mpart
00071       Integer                         :: inter  (2, ndim_3d, search%npart)
00072       Integer                         :: interl (2, ndim_3d, search%npart)
00073       Integer                         :: extent(ndim_3d)
00074       Integer                         :: idim(ndim_3d)
00075       Integer                         :: offset(ndim_3d)
00076       Integer                         :: addoffset(ndim_3d)
00077       Integer                         :: shift (ndim_3d)
00078       Integer, Allocatable            :: buffer (:)
00079 
00080 
00081 
00082       Type (Gridfunction), Pointer    :: field
00083       Integer                         :: n_vars_ret
00084 
00085 
00086 
00087       Integer                         :: defined, mask_id
00088       Logical                         :: src_mask_available, match
00089       Logical                         :: dst_mask_available
00090 
00091 
00092 
00093       Integer                         :: dir_index
00094 
00095 
00096 
00097       Integer                         :: ierrp (3)
00098 
00099 
00100 
00101 
00102 
00103 
00104 
00105 
00106 
00107 
00108 
00109 
00110 
00111 
00112 
00113 
00114 
00115 
00116 
00117 
00118    Character(len=len_cvs_string), save :: mycvs = 
00119        '$Id: psmile_search_donor_gridless.F90 2792 2010-12-01 17:35:05Z hanke $'
00120 
00121 
00122 
00123 
00124 
00125 #ifdef VERBOSE
00126       print 9990, trim(ch_id), search%msg_intersections%src_comp_id
00127 
00128       call psmile_flushstd
00129 #endif /* VERBOSE */
00130 
00131       ierror = 0
00132       n_vars_ret = 0
00133       npart = search%npart
00134       mpart = min(npart,maxpart)
00135 
00136 1000  continue
00137 
00138 #ifdef PRISM_ASSERTION
00139 
00140 
00141 
00142       if (        search%grid_type /= PRISM_Gridless .or. &
00143           Grids(grid_id)%grid_type /= PRISM_Gridless) then
00144           print *, trim(ch_id), ': search grid_type', &
00145                    search%grid_type, Grids(grid_id)%grid_type
00146           call psmile_assert ( __FILE__, __LINE__, &
00147                'Gridless Grids can coupled only with Gridless Grids')
00148       endif
00149 
00150       if (var_id < 1 .or. var_id > Number_of_Fields_allocated) then
00151          print *, 'var_id', var_id
00152          call psmile_assert (__FILE__, __LINE__, "Incorrect var_id")
00153       endif
00154 #endif
00155 
00156       if (npart > maxpart ) then
00157 
00158 #ifdef PRISM_ASSERTION
00159                   call psmile_assert (__FILE__, __LINE__, &
00160                                    "This part of the code does not work. I am sorry.")
00161 #endif
00162 #ifdef DEBUG
00163 
00164 
00165 #endif
00166 
00167 
00168 
00169 
00170 
00171 
00172 
00173 
00174 
00175 
00176 
00177 
00178 
00179 
00180       endif
00181 
00182 
00183 
00184 
00185 
00186 
00187       do ipart = 1, mpart
00188 
00189          inter(:,:,ipart) = search%msg_intersections%intersections(npart+ipart)%intersection
00190 
00191       end do
00192 
00193 
00194 
00195 
00196 
00197 
00198 
00199 
00200       shift = 0
00201 
00202       interl = inter
00203 
00204 
00205 
00206 
00207 
00208 
00209 
00210 
00211 
00212 
00213       if (Associated (Grids(grid_id)%partition)) then
00214 
00215          nbr_blocks = size(Grids(grid_id)%partition(:,1))
00216 
00217          if ( nbr_blocks > 1 ) then
00218 
00219             do n = 1, npart
00220 
00221                
00222 
00223                do nb = 1, nbr_blocks
00224 
00225                   if ( inter (2,1,n) >= Grids(grid_id)%partition(nb,1) + 1 .and. &
00226                        inter (1,1,n) <= Grids(grid_id)%partition(nb,1) +         &
00227                                         Grids(grid_id)%extent   (nb,1)     .and. &
00228                        inter (2,2,n) >= Grids(grid_id)%partition(nb,2) + 1 .and. &
00229                        inter (1,2,n) <= Grids(grid_id)%partition(nb,2) +         &
00230                                         Grids(grid_id)%extent   (nb,2)     .and. &
00231                        inter (2,3,n) >= Grids(grid_id)%partition(nb,3) + 1 .and. &
00232                        inter (1,3,n) <= Grids(grid_id)%partition(nb,3) +         &
00233                                         Grids(grid_id)%extent   (nb,3) ) exit
00234                enddo
00235 
00236                
00237 
00238                if ( nb <= nbr_blocks ) then
00239 #ifdef DEBUG
00240                   print *, ' Transform block ', nb, Grids(grid_id)%partition(nb,1)+1, &
00241                                          ' - ',     Grids(grid_id)%partition(nb,1)+ &
00242                                                     Grids(grid_id)%extent   (nb,1)
00243                   print *, ' global inter    ', nb, inter(:, :, n)
00244 #endif
00245                   offset    = 0
00246                   addoffset = 0
00247 
00248                   do i = 1, ndim_3d
00249                      idim(i) = Grids(grid_id)%grid_shape(2,i) - &
00250                                Grids(grid_id)%grid_shape(1,i)+1
00251                      do nl = 1, nb-1
00252                         addoffset(i) = Grids(grid_id)%extent(nl,i) + offset(i)
00253                         if ( addoffset(i) < idim(i) ) &
00254                              offset(i) = addoffset(i)
00255                      enddo
00256                   enddo
00257 
00258 
00259                   do i = 1, ndim_3d
00260                      extent(i) = interl(2,i,n) - interl(1,i,n)
00261                      interl(1, i, n) = interl(1, i, n) + offset(i) &
00262                                      - Grids(grid_id)%partition(nb,i) &
00263                                      + Grids(grid_id)%grid_shape(1,i) - 1
00264                      interl(2, i, n) = interl(1, i, n) + extent(i)
00265                   end do
00266 #ifdef DEBUG
00267                   print *, ' part offset ', nb, offset(:)
00268                   print *, ' local inter ', nb, interl(:, :, n)
00269 #endif
00270                else
00271 
00272                   ierror = PRISM_Error_Internal
00273                   ierrp (1) = nb
00274                   ierrp (2) = nbr_blocks
00275                   call psmile_error (ierror, "No block found", ierrp, 2, &
00276                        __FILE__, __LINE__)
00277                   return
00278 
00279                endif
00280 
00281             enddo
00282 
00283          else
00284 
00285             do n = 1, npart
00286                do i = 1, ndim_3d
00287                   interl (1:2, i, n) = &
00288                   interl (1:2, i, n) - Grids(grid_id)%partition (1,i) &
00289                                      + Grids(grid_id)%grid_shape(1,i) - 1
00290                end do
00291             end do
00292 
00293          endif
00294 
00295       endif
00296 
00297 
00298 
00299       field => Fields(var_id)
00300 
00301       mask_id = field%mask_id
00302       src_mask_available = mask_id /= PRISM_UNDEFINED
00303       dst_mask_available = Associated (search%search_mask)
00304 
00305       if (.not. dst_mask_available) then
00306 
00307 
00308 
00309 
00310 
00311          if (src_mask_available) then
00312             call psmile_is_mask_defined (Masks(mask_id)%mask_array, &
00313                                          Masks(mask_id)%mask_shape, &
00314                                          interl, npart, defined, ierror)
00315             if (ierror > 0) return
00316 
00317             if (defined /= 2) then
00318                ierrp (1) = grid_id
00319                ierrp (2) = search%msg_intersections%src_comp_id
00320                ierrp (3) = mask_id
00321                ierror = PRISM_Error_Mask
00322 
00323                call psmile_error ( ierror,                     &
00324                        "masks for gridless grids don't match", &
00325                        ierrp, 3, __FILE__, __LINE__ )
00326                return
00327             endif
00328          endif
00329 
00330 
00331 
00332 
00333 
00334          call psmile_locations_direct (interl, inter, search, method_id,  &
00335                                        dir_index, ierror)
00336          if (ierror > 0) return
00337       else 
00338 
00339 
00340 
00341          if (src_mask_available) then
00342                do ipart = 1, npart
00343                call psmile_do_masks_match (Masks(mask_id)%mask_array, &
00344                                            Masks(mask_id)%mask_shape, &
00345                                     search%search_mask(ipart)%vector, &
00346                                     search%shape(:,:,ipart),          &
00347                              inter(:, :, ipart), 1, match, ierror)
00348                if (ierror > 0) return
00349 
00350                if (.not. match) exit
00351                end do 
00352          else
00353             match = .false.
00354          endif
00355 
00356          if (.not. match) then
00357             ierrp (1) = grid_id
00358             ierrp (2) = search%msg_intersections%src_comp_id
00359             ierrp (3) = mask_id
00360             ierror = PRISM_Error_Mask
00361 
00362             call psmile_error ( ierror,                     &
00363                     "masks for gridless grids don't match", &
00364                     ierrp, 3, __FILE__, __LINE__ )
00365             return
00366          endif
00367 
00368 
00369 
00370          call psmile_locations_3d_mask (search, inter, shift, method_id, &
00371                                         dir_index, ierror)
00372          if (ierror > 0) return
00373 
00374       endif
00375 
00376 
00377 
00378       call psmile_store_send_info (var_id, search%msg_intersections%transient_out_id, &
00379                                    dir_index, cpl_index, PRISM_UNDEFINED, ierror)
00380       if (ierror > 0) return
00381 
00382 
00383 
00384       call psmile_return_locations_3d (search%msg_intersections,     &
00385                                        search%sender, method_id,     &
00386                                        dir_index, cpl_index, n_vars, &
00387                                        n_vars_ret, ierror)
00388       if (ierror > 0) return
00389 
00390 
00391 
00392       if (n_vars_ret < n_vars) then
00393 
00394          call psmile_get_next_field (comp_info, search, field_list, &
00395                      n_vars, n_vars_ret, var_id, ierror)
00396          if (ierror > 0) return
00397 
00398          method_id = Fields(var_id)%method_id
00399 
00400 
00401          go to 1000
00402       endif
00403 
00404 
00405 
00406 #ifdef VERBOSE
00407       print 9980, trim(ch_id), search%msg_intersections%src_comp_id, ierror
00408 
00409       call psmile_flushstd
00410 #endif /* VERBOSE */
00411 
00412 
00413 
00414 #ifdef VERBOSE
00415 
00416 9990 format (1x, a, ': psmile_search_donor_gridless: comp_id =', i3)
00417 9980 format (1x, a, ': psmile_search_donor_gridless: comp_id =', i3, &
00418                     '; eof ierror =', i4)
00419 
00420 #endif /* VERBOSE */
00421 
00422       end subroutine PSMILe_Search_donor_gridless