00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_search_donor_irreg2_dble (comp_info, &
00012 found, locations, len, search, &
00013 field_list, n_vars, &
00014 grid_id, method_id, var_id, tol, ierror)
00015
00016
00017
00018 use PRISM_constants
00019 use psmile_grid, only : common_grid_range
00020 use PSMILe, dummy_interface => PSMILe_Search_donor_irreg2_dble
00021
00022 implicit none
00023
00024
00025
00026 Type (Enddef_comp), Intent (In) :: comp_info
00027
00028
00029
00030
00031 Type (Enddef_search) :: search
00032
00033
00034
00035 Integer, Intent (In) :: len (search%npart, 2)
00036
00037
00038
00039
00040
00041
00042 Integer, Intent (In) :: grid_id
00043
00044
00045
00046 Integer, Intent (InOut) :: method_id
00047
00048
00049
00050 Integer, Intent (InOut) :: var_id
00051
00052
00053
00054 Integer, Intent (In) :: n_vars
00055
00056
00057
00058 Integer, Intent (In) :: field_list (nd_field_list, n_vars)
00059
00060
00061
00062 Double Precision, Intent (In) :: tol
00063
00064
00065
00066
00067
00068
00069 Type (integer_vector) :: found (search%npart, 2)
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080 Type (integer_vector) :: locations (search%npart, 2)
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090 integer, Intent (Out) :: ierror
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102 Integer, Parameter :: indl = 1
00103 Integer, Parameter :: indz = 2
00104 Integer, parameter :: nc_reg = 2
00105
00106 Logical, Parameter :: virtual_cell_required = .false.
00107
00108 Integer :: ipart, npart ,ndim
00109
00110
00111
00112 Integer :: i, j, ii, jj, j1, j2
00113 Type (integer_vector) :: found_3d (search%npart)
00114 Type (integer_vector) :: locations_3d (search%npart)
00115 Type (integer_vector) :: virtual_cell (search%npart)
00116 Type (integer_vector) :: virtual_3d (search%npart)
00117
00118
00119
00120 Type (dble_vector) :: array (ndim_3d, search%npart)
00121
00122 Integer :: shape_2d (2, ndim_3d, search%npart)
00123 Integer :: range_2d (2, ndim_3d, search%npart)
00124 Integer :: control_2d (2, ndim_3d, search%npart)
00125 Integer :: save_2d (2, ndim_3d, search%npart)
00126
00127 Integer :: tmp_range (2, ndim_3d, search%npart)
00128
00129 Integer :: shape_1d (2, ndim_3d, search%npart)
00130 Integer :: range_1d (2, ndim_3d, search%npart)
00131 Integer :: control_1d (2, ndim_3d, search%npart)
00132 Integer :: save_1d (2, ndim_3d, search%npart)
00133
00134 Integer :: dest
00135 Type (enddef_msg_locations) :: msg_locations
00136 Integer :: msgloc (msgloc_size)
00137
00138
00139
00140 Type (Corner_Block), Pointer :: corner_pointer
00141 Integer :: n_vars_ret
00142
00143
00144
00145 Integer, Parameter :: lev = 1
00146 Integer :: dim
00147 Integer :: method_type, old1, old2
00148 Integer :: cpl_index, dir_index
00149 Integer :: len_cpl (search%npart)
00150 Integer :: m_levdim (ndim_3d)
00151 Type (Coords_Block),Pointer :: coords_pointer
00152
00153 Logical :: changed
00154 Logical :: msk_required
00155 Type (integer_vector) :: locations_save (search%npart, 2)
00156 Type (integer_vector) :: found_save (search%npart, 2)
00157
00158 Integer :: lenl (search%npart, 2)
00159 Integer :: lenc (ndim_3d)
00160 Type (Enddef_mg_double) :: m_arrays
00161
00162
00163
00164 Integer :: nlev
00165
00166
00167
00168 Double Precision, Save :: period (ndim_2d)
00169
00170
00171
00172
00173 Integer, parameter :: nerrp = 2
00174 Integer :: ierrp (nerrp)
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199 Character(len=len_cvs_string), save :: mycvs =
00200 '$Id: psmile_search_donor_irreg2_dble.F90 2804 2010-12-07 10:07:10Z hanke $'
00201
00202
00203
00204
00205
00206 #ifdef VERBOSE
00207 print 9990, trim(ch_id), Grids(grid_id)%comp_id, search%sender
00208
00209 call psmile_flushstd
00210 #endif /* VERBOSE */
00211
00212 period (1:ndim_2d) = common_grid_range(2,1:ndim_2d) - &
00213 common_grid_range(1,1:ndim_2d)
00214
00215 n_vars_ret = 0
00216 npart = search%npart
00217 lenl = len
00218 old2 = PSMILe_Undef
00219
00220 corner_pointer => Grids(grid_id)%corner_pointer
00221
00222 #ifdef PRISM_ASSERTION
00223 if (method_id > Number_of_Methods_allocated .or. &
00224 method_id < 1) then
00225 print *, trim(ch_id), "PSMILe_Search_donor_irreg2_dble: method_id =", &
00226 method_id
00227 call psmile_assert (__FILE__, __LINE__, &
00228 "Invalid Method id")
00229 endif
00230 #endif
00231
00232
00233
00234
00235
00236
00237
00238 if (search%grid_type == PRISM_Reglonlatvrt) then
00239
00240
00241
00242 call psmile_srch_coords_reg_21d_dble (search, &
00243 array, shape_2d, range_2d, shape_1d, range_1d, &
00244 ierror)
00245
00246 else
00247
00248
00249
00250 do ipart = 1, npart
00251 array(1,ipart)%vector => search%search_dble(1,ipart)%vector
00252 array(2,ipart)%vector => search%search_dble(2,ipart)%vector
00253 array(3,ipart)%vector => search%search_dble(3,ipart)%vector
00254 end do
00255
00256
00257
00258 if (search%grid_type == PRISM_Irrlonlat_regvrt) then
00259
00260 shape_2d(:, 1:ndim_2d, 1:npart) = search%shape(:, 1:ndim_2d, 1:npart)
00261 shape_2d(:, ndim_3d, 1:npart) = 1
00262
00263 range_2d(:, 1:ndim_2d, 1:npart) = search%range(:, 1:ndim_2d, 1:npart)
00264 range_2d(:, ndim_3d, 1:npart) = 1
00265
00266 shape_1d(:, 1, 1:npart) = search%range(:, ndim_3d, 1:npart)
00267 shape_1d(:, 2:ndim_3d, 1:npart) = 1
00268
00269 range_1d = shape_1d
00270
00271 else if (search%grid_type == PRISM_Gaussreduced_regvrt) then
00272
00273 shape_2d(:, 1:ndim_2d, 1:npart) = search%shape(:, 1:ndim_2d, 1:npart)
00274 shape_2d(:, ndim_3d, 1:npart) = 1
00275
00276 range_2d(:, 1:ndim_2d, 1:npart) = search%range(:, 1:ndim_2d, 1:npart)
00277 range_2d(:, ndim_3d, 1:npart) = 1
00278
00279 shape_1d(:, 1, 1:npart) = search%range(:, ndim_3d, 1:npart)
00280 shape_1d(:, 2:ndim_3d, 1:npart) = 1
00281
00282 range_1d = shape_1d
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293 else
00294
00295 shape_2d(:, :, 1:npart) = search%shape(:, :, 1:npart)
00296 range_2d(:, :, 1:npart) = search%range(:, :, 1:npart)
00297
00298 shape_1d(:, :, 1:npart) = search%shape(:, :, 1:npart)
00299 range_1d(:, :, 1:npart) = search%range(:, :, 1:npart)
00300
00301 endif
00302
00303 endif
00304
00305 control_2d = range_2d
00306 control_1d = range_1d
00307
00308 save_2d = control_2d
00309 save_1d = control_1d
00310
00311
00312 nlev = Grids(grid_id)%nlev
00313
00314
00315
00316
00317
00318
00319
00320
00321 changed = .false.
00322
00323
00324
00325
00326
00327
00328
00329 if (n_vars > 0) then
00330 do ipart = 1, search%npart
00331 Allocate (locations_save(ipart,indl)%vector(1:len(ipart,indl)*ndim_2d), &
00332 locations_save(ipart,indz)%vector(1:len(ipart,indz)), &
00333 found_save(ipart,indl)%vector(1:len(ipart,indl)), &
00334 found_save(ipart,indz)%vector(1:len(ipart,indz)), STAT = ierror)
00335 if ( ierror > 0 ) then
00336 ierrp (1) = ierror
00337 ierrp (2) = 3*len(ipart,indl)+2*len(ipart,indz)
00338 ierror = PRISM_Error_Alloc
00339 call psmile_error ( ierror, 'locations_save and found_save', &
00340 ierrp, 2, __FILE__, __LINE__ )
00341 return
00342 endif
00343 end do
00344
00345 do ipart = 1, search%npart
00346
00347 locations_save(ipart,indl)%vector(1:len(ipart,indl)*ndim_2d) = &
00348 locations(ipart,indl)%vector(1:len(ipart,indl)*ndim_2d)
00349 locations_save(ipart,indz)%vector(1:len(ipart,indz)) = &
00350 locations(ipart,indz)%vector(1:len(ipart,indz))
00351
00352 found_save(ipart,indl)%vector(1:len(ipart,indl)) = &
00353 found(ipart,indl)%vector(1:len(ipart,indl))
00354 found_save(ipart,indz)%vector(1:len(ipart,indz)) = &
00355 found(ipart,indz)%vector(1:len(ipart,indz))
00356
00357 enddo
00358 endif
00359
00360
00361
00362 1000 continue
00363
00364 method_type = Methods(method_id)%method_type
00365 coords_pointer => Methods(method_id)%coords_pointer
00366
00367
00368
00369
00370 control_2d = save_2d
00371 control_1d = save_1d
00372
00373 if (method_type == PSMILe_PointMethod) then
00374
00375 if (coords_pointer%coords_datatype /= MPI_DOUBLE_PRECISION) then
00376 ierror = PRISM_Error_Internal
00377
00378 call psmile_error ( ierror, 'Different datatypes in Grid and Method are currently not supported', &
00379 ierrp, 0, __FILE__, __LINE__ )
00380 endif
00381
00382 #ifdef PRISM_ASSERTION
00383 if (.not. Associated(coords_pointer%coords_dble(1)%vector) ) then
00384 call psmile_assert (__FILE__, __LINE__, &
00385 "x coordinates of method are not defined")
00386 endif
00387 #endif
00388
00389
00390
00391
00392
00393 if ( old2 /= search%msg_intersections%requires_conserv_remap .and. changed ) then
00394
00395
00396
00397
00398
00399
00400 if ( search%grid_type == PRISM_Gaussreduced_regvrt ) then
00401
00402 if ( old2 == PSMILe_conserv2D ) then
00403 range_2d(:, 1:ndim_2d, 1:npart) = search%range(:, 1:ndim_2d, 1:npart)
00404
00405
00406 else
00407 ierrp (1) = search%grid_type
00408 ierror = PRISM_Error_Internal
00409 call psmile_error ( ierror, 'unsupported interpolation type for search grid', &
00410 ierrp, 1, __FILE__, __LINE__ )
00411 endif
00412
00413 else
00414
00415 do i = 1, ndim_3d
00416 do ipart = 1, search%npart
00417 if ( old2 == PSMILe_conserv3D ) then
00418 search%range(2,i,ipart) = search%range(2,i,ipart)-1
00419 shape_2d (2,i,ipart) = shape_2d (2,i,ipart)-1
00420 range_2d (2,i,ipart) = range_2d (2,i,ipart)-1
00421 range_1d (2,i,ipart) = range_2d (2,i,ipart)-1
00422 control_2d (2,i,ipart) = range_2d (2,i,ipart)
00423 control_1d (2,i,ipart) = range_1d (2,i,ipart)
00424 else if ( old2 == PSMILe_conserv2D .and. i < ndim_3d ) then
00425 search%range(2,i,ipart) = search%range(2,i,ipart)-1
00426 shape_2d (2,i,ipart) = shape_2d (2,i,ipart)-1
00427 range_2d (2,i,ipart) = range_2d (2,i,ipart)-1
00428 control_2d (2,i,ipart) = range_2d (2,i,ipart)
00429 else if ( old2 == PSMILe_conserv2D .and. i == ndim_3d ) then
00430 control_2d (2,i,ipart) = range_2d (2,i,ipart)
00431 control_1d (2,i,ipart) = range_1d (2,i,ipart)
00432 endif
00433 enddo
00434 enddo
00435
00436 endif
00437
00438 if ( search%grid_type == PRISM_Reglonlatvrt ) then
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450 do ipart = 1, search%npart
00451 j1 = shape_2d(2,1,ipart)-shape_2d(1,1,ipart) + 1
00452 do j = 1, shape_2d(2,2,ipart)-shape_2d(1,2,ipart) + 1
00453 array(1,ipart)%vector((j-1)*j1+1:j*j1) = &
00454 search%search_dble(1,ipart)%vector(1:j1)
00455 array(2,ipart)%vector((j-1)*j1+1:j*j1) = &
00456 search%search_dble(2,ipart)%vector(j)
00457 end do
00458 enddo
00459
00460 endif
00461
00462 do ipart = 1, search%npart
00463 lenl(ipart,indl) = ( search%range (2,1,ipart) - search%range (1,1,ipart) + 1 ) * &
00464 ( search%range (2,2,ipart) - search%range (1,2,ipart) + 1 )
00465 enddo
00466
00467 if ( old2 == PSMILe_conserv3D ) lenl(:,indz) = lenl (:,indz)-1
00468
00469 do i = 1, indz
00470 do ipart = 1, search%npart
00471
00472
00473
00474 Deallocate(found(ipart,i)%vector, locations(ipart,i)%vector, STAT = ierror)
00475
00476 if ( ierror > 0 ) then
00477 ierrp (1) = ierror
00478 ierrp (2) = m_levdim(i)
00479 ierror = PRISM_Error_Dealloc
00480 call psmile_error ( ierror, 'locations and found', &
00481 ierrp, 2, __FILE__, __LINE__ )
00482 return
00483 endif
00484
00485
00486
00487 if ( i == indl ) then
00488 Allocate (found (ipart,i)%vector(1:lenl(ipart,i)), &
00489 locations(ipart,i)%vector(1:lenl(ipart,i)*ndim_2d), STAT = ierror)
00490 else
00491 Allocate (found (ipart,i)%vector(1:lenl(ipart,i)), &
00492 locations(ipart,i)%vector(1:lenl(ipart,i)), STAT = ierror)
00493 endif
00494
00495 if ( ierror > 0 ) then
00496 ierrp (1) = ierror
00497 ierrp (2) = m_levdim(i)
00498 ierror = PRISM_Error_Alloc
00499 call psmile_error ( ierror, 'locations and found', &
00500 ierrp, 2, __FILE__, __LINE__ )
00501 return
00502 endif
00503
00504
00505
00506 if ( i == indl ) then
00507 j1 = 1
00508 j2 = 1
00509 do jj = search%range (1,2,ipart), search%range (2,2,ipart)
00510 do ii = search%range (1,1,ipart), search%range (2,1,ipart)
00511 locations(ipart,i)%vector(j1*2-1) = locations_save(ipart,i)%vector(j2*2-1)
00512 locations(ipart,i)%vector(j1*2) = locations_save(ipart,i)%vector(j2*2)
00513 found (ipart,i)%vector(j1) = found_save (ipart,i)%vector(j2)
00514 j1 = j1 + 1
00515 j2 = j2 + 1
00516 if ( ii == search%range (2,1,ipart) ) j2 = j2 + 1
00517 enddo
00518 enddo
00519 else
00520 locations(ipart,i)%vector = locations_save(ipart,i)%vector
00521 found (ipart,i)%vector = found_save (ipart,i)%vector
00522 endif
00523
00524 enddo
00525 enddo
00526
00527 else if (changed) then
00528
00529 do ipart = 1, search%npart
00530 found(ipart,indl)%vector(1:len(ipart,indl)) = &
00531 abs (found(ipart,indl)%vector(1:len(ipart,indl)))
00532 found(ipart,indz)%vector(1:len(ipart,indz)) = &
00533 abs (found(ipart,indz)%vector(1:len(ipart,indz)))
00534 enddo
00535
00536 do ipart = 1, search%npart
00537 locations (ipart,indl)%vector(1:len(ipart,indl)*ndim_2d) = &
00538 locations_save(ipart,indl)%vector(1:len(ipart,indl)*ndim_2d)
00539 locations (ipart,indz)%vector(1:len(ipart,indz)) = &
00540 locations_save(ipart,indz)%vector(1:len(ipart,indz))
00541 enddo
00542
00543 endif
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554 m_levdim (1:ndim_3d) = Grids(grid_id)%grid_shape(2,1:ndim_3d) - &
00555 Grids(grid_id)%grid_shape(1,1:ndim_3d) + 2
00556
00557
00558
00559 dim = m_levdim(1) * m_levdim (2)
00560
00561 do i = 1, ndim_2d
00562
00563 Allocate (m_arrays%chmin(i)%vector(dim), STAT = ierror)
00564 if ( ierror > 0 ) then
00565 ierrp (1) = ierror
00566 ierrp (2) = dim
00567 ierror = PRISM_Error_Alloc
00568 call psmile_error ( ierror, 'm_arrays%chmin(i)%vector', &
00569 ierrp, 2, __FILE__, __LINE__ )
00570 return
00571 endif
00572
00573 Allocate (m_arrays%chmax(i)%vector(dim), STAT = ierror)
00574 if ( ierror > 0 ) then
00575 ierrp (1) = ierror
00576 ierrp (2) = dim
00577 ierror = PRISM_Error_Alloc
00578 call psmile_error ( ierror, 'm_arrays%chmax(i)%vector', &
00579 ierrp, 2, __FILE__, __LINE__ )
00580 return
00581 endif
00582
00583 Allocate (m_arrays%midp(i)%vector(dim), STAT = ierror)
00584 if ( ierror > 0 ) then
00585 ierrp (1) = ierror
00586 ierrp (2) = dim
00587 ierror = PRISM_Error_Alloc
00588 call psmile_error ( ierror, 'm_arrays%midp(i)%vector', &
00589 ierrp, 2, __FILE__, __LINE__ )
00590 return
00591 endif
00592
00593 end do
00594
00595 i = ndim_3d
00596 Allocate (m_arrays%chmin(i)%vector(m_levdim(i)), STAT = ierror)
00597 if ( ierror > 0 ) then
00598 ierrp (1) = ierror
00599 ierrp (2) = m_levdim(i)
00600 ierror = PRISM_Error_Alloc
00601 call psmile_error ( ierror, 'm_arrays%chmin(3)%vector', &
00602 ierrp, 2, __FILE__, __LINE__ )
00603 return
00604 endif
00605
00606 Allocate (m_arrays%chmax(i)%vector(m_levdim(i)), STAT = ierror)
00607 if ( ierror > 0 ) then
00608 ierrp (1) = ierror
00609 ierrp (2) = m_levdim(i)
00610 ierror = PRISM_Error_Alloc
00611 call psmile_error ( ierror, 'm_arrays%chmax(3)%vector', &
00612 ierrp, 2, __FILE__, __LINE__ )
00613 return
00614 endif
00615
00616
00617
00618 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D .or. &
00619 search%msg_intersections%requires_conserv_remap == PSMILe_conserv3D ) then
00620
00621 lenc(:) = corner_pointer%corner_shape(2,:) - &
00622 corner_pointer%corner_shape(1,:) + 1
00623
00624 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D ) &
00625 ndim = ndim_2d
00626 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv3D ) &
00627 ndim = ndim_3d
00628
00629 if ( ndim == ndim_2d ) then
00630 i = ndim_3d
00631 call psmile_bbcells_1d_dble ( &
00632 coords_pointer%coords_dble(i)%vector, &
00633 coords_pointer%coords_shape(1:2, i), &
00634 Grids(grid_id)%grid_shape (1:2, i), &
00635 corner_pointer%corners_dble(i)%vector, &
00636 corner_pointer%corner_shape(1:2,i), nc_reg, &
00637 m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector, &
00638 m_levdim(i), ierror)
00639 if (ierror > 0) return
00640 endif
00641
00642 do ipart = 1, npart
00643 call psmile_mg_cells_2d_dble ( grid_id, search%grid_type, &
00644 found (ipart, indl)%vector, &
00645 locations(ipart, indl)%vector, range_2d(:, :, ipart), &
00646 array(1,ipart)%vector, array(2,ipart)%vector, &
00647 shape_2d(:, :, ipart), control_2d(:, :, ipart), &
00648 Grids(grid_id)%grid_shape(:,1:ndim_2d), ipart, &
00649 corner_pointer%corner_shape(:,1:ndim_2d), &
00650 corner_pointer%nbr_corners/nc_reg, &
00651 corner_pointer%corners_dble(1)%vector, &
00652 corner_pointer%corners_dble(2)%vector, &
00653 m_arrays%chmin(1)%vector, m_arrays%chmax(1)%vector, &
00654 m_arrays%chmin(2)%vector, m_arrays%chmax(2)%vector, &
00655 tol, ierror)
00656 if (ierror > 0) return
00657
00658
00659
00660
00661 if ( ndim == ndim_2d ) then
00662
00663 i = ndim_3d
00664 call psmile_mg_method_1d_dble ( comp_info, Grids(grid_id)%nlev, &
00665 found (ipart, indz)%vector, &
00666 locations(ipart, indz)%vector, range_1d(:, :, ipart), &
00667 search%search_dble(i, ipart)%vector, &
00668 shape_1d(:, :, ipart), control_1d(:, :, ipart), &
00669 method_id, &
00670 coords_pointer%coords_dble(i)%vector, &
00671 coords_pointer%coords_shape(:, i), &
00672 Grids(grid_id)%grid_shape (:, i), &
00673 Grids(grid_id)%cyclic(i), &
00674 m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector, &
00675 tol, ierror)
00676 if (ierror > 0) return
00677
00678 else if ( ndim == ndim_3d ) then
00679
00680 ierror = PSMILe_Undef
00681 call psmile_error ( ierror, 'Interpolation type not supported', &
00682 ierrp, 0, __FILE__, __LINE__ )
00683
00684 call psmile_mg_cells_1d_dble ( Grids(grid_id)%nlev, &
00685 found (ipart, indz)%vector, &
00686 locations(ipart, indz)%vector, range_2d(:, :, ipart), &
00687 search%grid_type, &
00688 search%search_dble(i, ipart)%vector, i, &
00689 shape_1d(:,:,ipart), control_1d(:, :, ipart), &
00690 Grids(grid_id)%grid_shape (:, i), &
00691 Grids(grid_id)%cyclic(i), &
00692 m_arrays%chmin(i)%vector, m_arrays%chmax(i)%vector, &
00693 tol, ierror)
00694 if (ierror > 0) return
00695
00696 endif
00697
00698 enddo
00699
00700 else
00701
00702
00703
00704 call psmile_mg_method_irreg2_dble (comp_info, &
00705 found, locations, search, &
00706 array, shape_2d, range_2d, control_2d, &
00707 shape_1d, range_1d, control_1d, &
00708 m_arrays, m_levdim, &
00709 grid_id, method_id, tol, ierror)
00710
00711 endif
00712
00713 changed = .true.
00714
00715
00716
00717 do i = 1, ndim_2d
00718
00719 Deallocate (m_arrays%midp (i)%vector)
00720 Deallocate (m_arrays%chmin(i)%vector)
00721 Deallocate (m_arrays%chmax(i)%vector)
00722 end do
00723
00724 i = ndim_3d
00725 Deallocate (m_arrays%chmin(i)%vector)
00726 Deallocate (m_arrays%chmax(i)%vector)
00727
00728 else if (method_type == PSMILe_VectorPointMethod) then
00729
00730 ierrp (1) = method_type
00731 ierror = PRISM_Error_Internal
00732 call psmile_error ( ierror, 'Vector Method is currently not supported', &
00733 ierrp, 1, __FILE__, __LINE__ )
00734
00735 else if (method_type == PSMILe_SubgridMethod) then
00736 ierrp (1) = method_type
00737 ierror = PRISM_Error_Internal
00738 call psmile_error ( ierror, 'Subgrid Method is currently not supported', &
00739 ierrp, 1, __FILE__, __LINE__ )
00740
00741 else
00742 ierrp (1) = method_type
00743 ierror = PRISM_Error_Internal
00744 call psmile_error ( ierror, 'unsupported method type', &
00745 ierrp, 1, __FILE__, __LINE__ )
00746
00747 endif
00748
00749
00750
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768 1500 continue
00769
00770
00771
00772 call psmile_compact_locations ( grid_id, search, ndim_2d, found, ierror )
00773
00774 if ( old2 == search%msg_intersections%requires_conserv_remap ) then
00775 control_2d = save_2d
00776 control_1d = save_1d
00777 endif
00778
00779 if (search%grid_type == PRISM_Irrlonlatvrt) then
00780
00781 control_2d (1, 1:ndim_3d, 1:npart) = &
00782 max (control_2d (1, 1:ndim_3d, 1:npart), &
00783 control_1d (1, 1:ndim_3d, 1:npart))
00784 control_2d (2, 1:ndim_3d, 1:npart) = &
00785 min (control_2d (2, 1:ndim_3d, 1:npart), &
00786 control_1d (2, 1:ndim_3d, 1:npart))
00787
00788 else
00789
00790 control_2d (:, ndim_3d, 1:npart) = &
00791 control_1d (:, 1, 1:npart)
00792
00793 endif
00794
00795 if ( search%grid_type == PRISM_Gaussreduced_regvrt ) then
00796 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D .or. &
00797 search%msg_intersections%requires_conserv_remap == PSMILe_conserv3D ) then
00798
00799
00800
00801 msk_required = .true.
00802 else
00803 msk_required = .false.
00804 endif
00805 else
00806 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D .or. &
00807 search%msg_intersections%requires_conserv_remap == PSMILe_conserv3D ) then
00808 control_2d (2,1:ndim,1:npart) = control_2d (2,1:ndim,1:npart) - 1
00809 msk_required = .true.
00810 else
00811 msk_required = .false.
00812 endif
00813 endif
00814
00815 if (search%grid_type == PRISM_Irrlonlatvrt .or. &
00816 Associated(search%search_mask)) then
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828 call psmile_mg_found_loc_to_3d (search, nlev, &
00829 PRISM_Irrlonlat_Regvrt, found, locations, lenl, &
00830 virtual_cell, virtual_cell_required, &
00831 found_3d, locations_3d, virtual_3d, ierror)
00832 if (ierror > 0) return
00833
00834
00835
00836 call psmile_locations_3d (found_3d, locations_3d, search%range, &
00837 control_2d, search, method_id, msk_required, &
00838 virtual_3d, virtual_cell_required, &
00839 dir_index, cpl_index, len_cpl, ierror)
00840 if (ierror > 0) return
00841
00842 do ipart = 1, npart
00843 Deallocate (locations_3d(ipart)%vector)
00844 Deallocate (found_3d (ipart)%vector)
00845 end do
00846
00847
00848
00849
00850
00851
00852
00853 else if (search%grid_type == PRISM_Gaussreduced_Regvrt) then
00854
00855 shape_2d (:, ndim_3d, 1:npart) = &
00856 shape_1d (:, 1, 1:npart)
00857
00858 tmp_range = search%range
00859
00860 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D ) then
00861
00862 else
00863 tmp_range = search%range
00864 endif
00865
00866
00867
00868 call psmile_locations_irreg2 (found, locations, tmp_range, &
00869 control_2d, search, method_id, msk_required, &
00870 virtual_cell, virtual_cell_required, &
00871 dir_index, cpl_index, len_cpl, ierror)
00872
00873 else
00874
00875 shape_2d (:, ndim_3d, 1:npart) = &
00876 shape_1d (:, 1, 1:npart)
00877
00878
00879
00880 call psmile_locations_irreg2 (found, locations, search%range, &
00881 control_2d, search, method_id, msk_required, &
00882 virtual_cell, virtual_cell_required, &
00883 dir_index, cpl_index, len_cpl, ierror)
00884
00885 endif
00886
00887 if (ierror > 0) return
00888
00889 if (cpl_index /= PRISM_UNDEFINED) then
00890
00891
00892 call psmile_info_trs_loc_irreg2_dble (comp_info, &
00893 array, shape_2d, control_2d, len_cpl, &
00894 var_id, Grids(grid_id)%grid_shape, &
00895 search, method_id, &
00896 cpl_index, ierror)
00897 if (ierror > 0) return
00898
00899 else
00900
00901 if ( search%msg_intersections%requires_conserv_remap == PSMILe_conserv2D ) then
00902
00903 dest = search%sender
00904
00905 call psmile_init_enddef_msg_locs (msg_locations)
00906
00907 msg_locations%requires_conserv_remap = search%msg_intersections%requires_conserv_remap
00908 msg_locations%src_rank = psmile_rank
00909 msg_locations%msg_len = 0
00910 msg_locations%transi_out_id = search%msg_intersections%transient_out_id
00911 msg_locations%transi_in_id = search%msg_intersections%transient_in_id
00912
00913 call psmile_pack_msg_locations (msg_locations, msgloc)
00914 #ifdef DEBUG
00915 print *, ' Bsend zero msg with tag ', &
00916 loctag+search%msg_intersections%relative_msg_tag, ' to ', dest
00917 #endif /* DEBUG */
00918 call psmile_bsend ( msgloc, msgloc_size, MPI_INTEGER, dest, &
00919 loctag+search%msg_intersections%relative_msg_tag, &
00920 comm_psmile, ierror )
00921
00922 endif
00923
00924 endif
00925
00926
00927
00928 call psmile_store_send_info (var_id, search%msg_intersections%transient_out_id, &
00929 dir_index, cpl_index, PRISM_UNDEFINED, ierror)
00930 if (ierror > 0) return
00931
00932
00933
00934 call psmile_return_locations_3d (search%msg_intersections, &
00935 search%sender, method_id, &
00936 dir_index, cpl_index, &
00937 n_vars, n_vars_ret, ierror)
00938 if (ierror > 0) return
00939
00940
00941
00942
00943 if (n_vars_ret < n_vars) then
00944
00945 old2 = search%msg_intersections%requires_conserv_remap
00946
00947 call psmile_get_next_field (comp_info, search, field_list, n_vars, &
00948 n_vars_ret, var_id, ierror)
00949
00950 if (ierror > 0) return
00951
00952 old1 = method_id
00953 method_id = Fields(var_id)%method_id
00954
00955
00956
00957
00958
00959
00960
00961 if (old1 == method_id .and. old2 == search%msg_intersections%requires_conserv_remap ) &
00962 go to 1500
00963 go to 1000
00964 endif
00965
00966
00967
00968 if (n_vars > 0) then
00969 do ipart = 1, search%npart
00970 Deallocate (locations_save(ipart,indl)%vector)
00971 Deallocate (locations_save(ipart,indz)%vector)
00972 enddo
00973 endif
00974
00975
00976
00977
00978
00979 if (search%grid_type == PRISM_Reglonlatvrt) then
00980 do ipart = 1, npart
00981 Deallocate (array(2,ipart)%vector)
00982 Deallocate (array(1,ipart)%vector)
00983 end do
00984 endif
00985
00986 #ifdef VERBOSE
00987 print 9980, trim(ch_id), Grids(grid_id)%comp_id, search%sender, ierror
00988
00989 call psmile_flushstd
00990 #endif /* VERBOSE */
00991
00992
00993
00994 9990 format (1x, a, ': psmile_search_donor_irreg2_dble: comp_id =', i3, &
00995 '; sender =', i4)
00996 9980 format (1x, a, ': psmile_search_donor_irreg2_dble: comp_id =', i3, &
00997 '; eof sender =', i3, ', ierror =', i4)
00998
00999 end subroutine PSMILe_Search_donor_irreg2_dble