00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_info_trs_locs_3d_dble (comp_info, &
00012 coords, shape, control, len_cpl, &
00013 var_id, grid_valid_shape, &
00014 search, method_id, &
00015 send_index, ierror)
00016
00017
00018
00019 use PRISM_constants
00020
00021 use PSMILe, dummy_interface => PSMILe_info_trs_locs_3d_dble
00022 #ifdef DEBUG_TRACE
00023 use psmile_debug_trace
00024 #endif
00025
00026 implicit none
00027
00028
00029
00030 Type (Enddef_comp), Intent (In) :: comp_info
00031
00032
00033
00034
00035 Type (Enddef_search), Intent (InOut) :: search
00036
00037
00038
00039 Type (dble_vector), Intent (In) :: coords (ndim_3d, search%npart)
00040
00041
00042
00043 Integer, Intent (In) :: shape (2, ndim_3d, search%npart)
00044
00045
00046
00047
00048
00049 Integer, Intent (In) :: control (2, ndim_3d, search%npart)
00050
00051
00052
00053 Integer, Intent (In) :: len_cpl (search%npart)
00054
00055
00056
00057 Integer, Intent (In) :: var_id
00058
00059
00060
00061 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00062
00063
00064
00065 Integer, Intent (In) :: method_id
00066
00067
00068
00069 Integer, Intent (InOut) :: send_index
00070
00071
00072
00073
00074
00075 integer, Intent (Out) :: ierror
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085 Integer, Parameter :: interpolation_mode = 3
00086 Integer, Parameter :: n_intmethods = 1
00087
00088
00089
00090
00091
00092 Type (Gridfunction), Pointer :: field
00093
00094
00095
00096 Type (Method), Pointer :: mp
00097
00098
00099
00100 Type (Grid), Pointer :: grid_info
00101 Integer :: grid_id
00102
00103 Integer :: ncpl
00104
00105
00106
00107 Integer :: i, n
00108
00109
00110
00111 Type (Send_information), Pointer :: send_info
00112 Type (Extra_search_info) :: extra_search
00113
00114 Integer :: nb_neighbors, nb_extra
00115
00116 Integer, Pointer :: dstijk (:, :)
00117 Integer, Pointer :: srcloc (:)
00118
00119 Integer, Allocatable :: neighbors_3d (:, :, :)
00120 Integer :: interpolation_type
00121 Integer :: interpolation_methods (n_intmethods)
00122 Logical :: interpolation_search (n_intmethods)
00123 Integer :: interpolation(ndim_3d)
00124
00125 Logical :: global_search
00126
00127
00128
00129
00130
00131
00132
00133 Integer :: id_src_mask, id_tgt_mask
00134 Integer :: mask_id
00135 Integer :: use_how (3)
00136 Logical :: src_mask_available, use_mask
00137
00138 Logical, target :: dummy_mask_array (1)
00139
00140 Integer :: dummy_mask_shape (2, ndim_3d)
00141 Integer, target :: dummy_mask (1)
00142
00143
00144 Integer :: mask_shape (2, ndim_3d)
00145 Logical, Pointer :: mask_array (:)
00146 Integer, Pointer :: src_mask (:)
00147
00148
00149
00150 Integer :: ibeg, iend, ipart, nprev
00151
00152
00153
00154 Integer :: nb_corners
00155
00156 Integer :: cpl_id
00157 Integer :: epio_id, trs_rank
00158 Integer :: local_src_size, src_size
00159 Integer :: id_trans_out
00160 Integer :: ip
00161 Logical :: allocate_src
00162
00163 Type (dble_vector) :: src_grid (ndim_3d)
00164 Integer, Allocatable :: neighbors (:, :)
00165
00166
00167
00168 Integer :: id_trans_in, dest_comp_id
00169 Logical :: use_target_grid
00170
00171 Type (dble_vector) :: tgt_grid (ndim_3d)
00172
00173
00174
00175 logical :: transformed
00176 Type (dble_vector) :: sinvec, cosvec
00177
00178
00179
00180 Integer, parameter :: nerrp = 3
00181 Integer :: ierrp (nerrp)
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205 Character(len=len_cvs_string), save :: mycvs =
00206 '$Id: psmile_info_trs_locs_3d_dble.F90 2927 2011-01-28 14:04:12Z hanke $'
00207
00208
00209
00210
00211
00212 data dummy_mask_shape /0, 1, 0, 1, 0, 1/
00213
00214
00215
00216
00217
00218 #ifdef VERBOSE
00219 print 9990, trim(ch_id), method_id, search%msg_intersections%first_tgt_method_id, &
00220 search%sender
00221
00222 call psmile_flushstd
00223 #endif /* VERBOSE */
00224
00225 ierror = 0
00226
00227 field => Fields(var_id)
00228 mp => Methods(method_id)
00229 send_info => mp%send_infos_coupler(send_index)
00230 grid_id = mp%grid_id
00231 grid_info => Grids(grid_id)
00232
00233 transformed = .false.
00234
00235 global_search = .false.
00236 interpolation_search (:) = .false.
00237
00238 use_how(1) = PSMILe_undef
00239 use_how(2) = PSMILe_undef
00240 use_how(3) = grid_id
00241 #ifdef PRISM_ASSERTION
00242
00243
00244
00245 if ( .not. associated(field%Taskout) ) then
00246 print *, 'var_id', var_id
00247 call psmile_assert (__FILE__, &
00248 __LINE__, "Taskout is not available for Field")
00249 end if
00250
00251 if (var_id < 1 .or. var_id > Number_of_Fields_allocated) then
00252 print *, 'var_id', var_id
00253 call psmile_assert (__FILE__, &
00254 __LINE__, "Incorrect var_id")
00255 endif
00256
00257 if ( mp%status == PSMILe_Status_free .or. &
00258 mp%status == PSMILe_Status_undefined ) then
00259 call psmile_assert (__FILE__, &
00260 __LINE__, "Incorrect method")
00261 endif
00262
00263 do ipart = 1, search%npart
00264 if (control (1,1,ipart) < search%range (1,1,ipart) .or. &
00265 search%range (2,1,ipart) < control (2,1,ipart) .or. &
00266 control (1,2,ipart) < search%range (1,2,ipart) .or. &
00267 search%range (2,2,ipart) < control (2,2,ipart) .or. &
00268 control (1,3,ipart) < search%range (1,3,ipart) .or. &
00269 search%range (2,3,ipart) < control (2,3,ipart)) then
00270 print *, ipart, control (:, :, ipart), search%range (:, :, ipart)
00271 call psmile_assert (__FILE__, &
00272 __LINE__, "control is out of range")
00273 endif
00274 end do
00275
00276 if (send_index < 1) then
00277 print *, send_index
00278 call psmile_assert (__FILE__, &
00279 __LINE__, "send_index is out of range")
00280 endif
00281
00282 if (send_info%nloc < 1) then
00283 print *, "ncpl", send_info%nloc
00284 call psmile_assert (__FILE__, &
00285 __LINE__, "ncpl > 0 expected")
00286 endif
00287
00288 if (send_info%nvec /= 1 ) then
00289 print *, "nvec", send_info%nvec
00290 call psmile_assert (__FILE__, &
00291 __LINE__, "nvec == 1 expected")
00292 endif
00293
00294 if (send_info%nparts /= 1 ) then
00295 print *, "nparts", send_info%nparts
00296 call psmile_assert (__FILE__, &
00297 __LINE__, "nparts == 1 expected")
00298 endif
00299
00300 if (send_info%nloc /= Sum(len_cpl)) then
00301 print *, "ncpl", send_info%nloc, Sum (len_cpl)
00302 call psmile_assert (__FILE__, &
00303 __LINE__, "ncpl == SUM(len_cpl) expected")
00304 endif
00305 #endif
00306
00307
00308
00309 ncpl = send_info%nloc
00310 dstijk => send_info%dstijk
00311
00312 srcloc => send_info%srclocs(1,1)%vector
00313
00314 #ifdef DEBUG_TRACE
00315
00316 do ictl = 1, ncpl
00317 if (dstijk (1, ictl) == ictl_ind (1) .and. &
00318 dstijk (2, ictl) == ictl_ind (2) .and. &
00319 dstijk (3, ictl) == ictl_ind (3)) exit
00320 end do
00321 #endif
00322
00323
00324
00325
00326
00327
00328
00329 mask_id = field%mask_id
00330 src_mask_available = mask_id /= PRISM_UNDEFINED
00331 use_mask = src_mask_available
00332
00333 if (src_mask_available) then
00334 mask_array => Masks(mask_id)%mask_array
00335
00336 mask_shape = Masks(mask_id)%mask_shape
00337 else
00338 mask_array => dummy_mask_array
00339
00340 mask_shape = dummy_mask_shape
00341 endif
00342
00343
00344
00345
00346
00347 id_trans_in = search%msg_intersections%transient_in_id
00348 id_trans_out = search%msg_intersections%transient_out_id
00349 dest_comp_id = all_comp_infos (search%msg_intersections%all_comp_infos_comp_idx)%global_comp_id
00350
00351 do i = 1, size (field%Taskout)
00352 if (field%Taskout(i)%origin_type == PSMILe_comp .and. &
00353 field%Taskout(i)%remote_comp_id == dest_comp_id .and. &
00354 field%Taskout(i)%remote_transi_id == id_trans_in .and. &
00355 field%Taskout(i)%global_transi_id == id_trans_out) exit
00356 end do
00357
00358 if (i > size (field%Taskout)) then
00359 #ifdef DEBUG
00360 print *, 'psmile_info_trs_locs_3d_dble: method_id ', method_id
00361 print *, 'search for Taskout id ', id_trans_out, ' and Taskin Id ', id_trans_in, &
00362 ' for dest component', dest_comp_id
00363
00364 do i = 1, size (field%Taskout)
00365 print 9960, i, &
00366 field%Taskout(i)%origin_type, &
00367 field%Taskout(i)%remote_comp_id, &
00368 field%Taskout(i)%remote_transi_id, &
00369 field%Taskout(i)%global_transi_id
00370 end do
00371 #endif /* DEBUG */
00372
00373 ierrp (1) = id_trans_in
00374 ierrp (2) = id_trans_out
00375 ierror = PRISM_Error_Internal
00376
00377 call psmile_error ( ierror, 'cannot find Taskout index for transient ids', &
00378 ierrp, 2, __FILE__, __LINE__ )
00379 endif
00380
00381 interpolation_type = field%Taskout(i)%interp%interp_type
00382 interpolation_methods (1) = field%Taskout(i)%interp%interp_meth (1)
00383
00384
00385
00386
00387
00388
00389 #ifdef DEBUG
00390 print *, "interpolation type ", interpolation_type
00391 print *, "interpolation method", interpolation_methods
00392 #endif
00393
00394 if (interpolation_type /= PSMILe_3D) then
00395 ierrp (1) = interpolation_type
00396
00397 ierror = PRISM_Error_Interp_type
00398 call psmile_error ( ierror, 'PRISM_Irrlonlatvrt (which requires 3d-interpolation)', &
00399 ierrp, 1, __FILE__, __LINE__ )
00400 return
00401 endif
00402
00403 if (interpolation_methods (1) == PSMILe_none) then
00404 ierror = PRISM_Error_Interp_type
00405 call psmile_error ( ierror, 'Interpolation is PSMILE_none', &
00406 ierrp, 0, __FILE__, __LINE__ )
00407 return
00408 endif
00409
00410
00411
00412 select case (interpolation_methods (1))
00413
00414
00415
00416 case (PSMILe_nnghbr3D)
00417 nb_neighbors = field%Taskout(i)%interp%arg2(1)
00418 nb_extra = nb_neighbors
00419 use_mask = use_mask .and. &
00420 (field%Taskout(i)%interp%arg5(1) /= 1)
00421 #ifdef DEBUG
00422 print *, "Number of nneigbours", nb_neighbors
00423 #endif
00424
00425 case (PSMILe_trilinear)
00426 nb_neighbors = 8
00427 use_how(1) = field%Taskout(i)%interp%arg4(1)
00428 interpolation_search (1) = &
00429 field%Taskout(i)%interp%arg3 (1) == PSMILE_global
00430
00431 nb_extra = 1
00432
00433 case (PSMILe_none)
00434
00435 case (PSMILe_undef)
00436 ierrp (1) = 1
00437 ierror = PRISM_Error_Interp_type
00438
00439 call psmile_error ( ierror, &
00440 'undefined interpolation method (PSMILE_UNDEF) found', &
00441 ierrp, 1, __FILE__, __LINE__ )
00442 return
00443
00444 case default
00445 ierrp (1) = interpolation_methods (1)
00446 ierror = PRISM_Error_Internal
00447
00448 call psmile_error ( ierror, &
00449 'unsupported interpolation method for type of grid PRISM_Irrlonlatvrt', &
00450 ierrp, 1, __FILE__, __LINE__ )
00451 return
00452
00453 end select
00454
00455 if (nb_neighbors < 1) then
00456 ierror = PRISM_Error_Interp_type
00457 ierrp (1) = nb_neighbors
00458
00459 call psmile_error ( ierror, &
00460 'Number of interpolation bases is less than 1', &
00461 ierrp, 1, __FILE__, __LINE__ )
00462 return
00463 endif
00464
00465 if (comp_info%size > 1) then
00466 global_search = interpolation_search (1)
00467
00468 #ifdef DEBUG
00469 print *, "interpolation search", interpolation_search (:)
00470 print *, "global search", global_search
00471 #endif
00472 else
00473 interpolation_search = .false.
00474 endif
00475
00476
00477
00478 Allocate (neighbors_3d(ndim_3d, ncpl, nb_neighbors), STAT = ierror)
00479 if ( ierror > 0 ) then
00480 ierrp (1) = ierror
00481 ierrp (2) = ndim_3d * ncpl * nb_neighbors
00482
00483 ierror = PRISM_Error_Alloc
00484 call psmile_error ( ierror, 'neighbors_3d', &
00485 ierrp, 2, __FILE__, __LINE__ )
00486 return
00487 endif
00488
00489
00490
00491
00492
00493 extra_search%global_search = global_search
00494
00495 use_target_grid = search%npart == 1 .and. &
00496 ncpl == (shape(2,1,1)-shape(1,1,1)+1)* &
00497 (shape(2,2,1)-shape(1,2,1)+1)* &
00498 (shape(2,3,1)-shape(1,3,1)+1)
00499
00500
00501
00502 if (use_target_grid) then
00503
00504
00505
00506 tgt_grid(1)%vector => coords(1,1)%vector
00507 tgt_grid(2)%vector => coords(2,1)%vector
00508 tgt_grid(3)%vector => coords(3,1)%vector
00509
00510 else
00511
00512
00513
00514 do i = 1, ndim_3d
00515 Allocate (tgt_grid(i)%vector(ncpl), stat = ierror)
00516 if ( ierror > 0 ) then
00517 ierrp (1) = ierror
00518 ierrp (2) = ncpl
00519
00520 ierror = PRISM_Error_Alloc
00521 call psmile_error ( ierror, 'tgt_grid(i)%vector', &
00522 ierrp, 2, __FILE__, __LINE__ )
00523 return
00524 endif
00525 end do
00526
00527
00528
00529 nprev = 0
00530 do ipart = 1, search%npart
00531 if (len_cpl (ipart) > 0) then
00532 ibeg = nprev + 1
00533 iend = nprev + len_cpl (ipart)
00534
00535 do i = 1, ndim_3d
00536 call psmile_extract_indices_3d_dble ( &
00537 coords(i,ipart)%vector, shape(:,:, ipart), &
00538 dstijk (:, ibeg:iend), len_cpl(ipart), &
00539 tgt_grid(i)%vector(ibeg:iend), &
00540 ierror)
00541 if (ierror > 0) return
00542 end do
00543
00544 nprev = nprev + len_cpl (ipart)
00545 endif
00546 end do
00547 endif
00548
00549 call psmile_neigh_extra_search_init (search, grid_id, extra_search, ierror)
00550 if (ierror > 0) return
00551
00552
00553
00554 select case (interpolation_methods (1))
00555
00556
00557
00558 case (PSMILe_nnghbr3D)
00559
00560
00561
00562
00563 call psmile_info_trf_coords_3d_dble ( &
00564 mp%coords_pointer%coords_dble(1)%vector, &
00565 mp%coords_pointer%coords_dble(2)%vector, &
00566 mp%coords_pointer%coords_dble(3)%vector, &
00567 mp%coords_pointer%coords_shape, &
00568 grid_valid_shape, &
00569 sinvec, cosvec, ierror)
00570 if (ierror > 0) return
00571
00572 transformed = .true.
00573
00574
00575
00576 nprev = 0
00577
00578 do ipart = 1, search%npart
00579 if (len_cpl (ipart) > 0) then
00580
00581 call psmile_neigh_nearest_3d_dble (grid_id, &
00582 tgt_grid(1)%vector, tgt_grid(2)%vector, &
00583 tgt_grid(3)%vector, &
00584 mp%coords_pointer%coords_dble(1)%vector, &
00585 mp%coords_pointer%coords_dble(2)%vector, &
00586 mp%coords_pointer%coords_dble(3)%vector, &
00587 mp%coords_pointer%coords_shape, &
00588 mask_array, mask_shape, use_mask, &
00589 sinvec%vector, cosvec%vector, grid_valid_shape, &
00590 srcloc, ncpl, nprev, len_cpl (ipart), &
00591 neighbors_3d, nb_neighbors, &
00592 extra_search, ierror)
00593
00594 if (ierror > 0) return
00595
00596 nprev = nprev + len_cpl (ipart)
00597 end if
00598 end do
00599
00600
00601
00602 case (PSMILe_trilinear)
00603
00604 call psmile_neigh_trili_3d (grid_valid_shape, &
00605 interpolation_mode, grid_info%cyclic, &
00606 srcloc, ncpl, &
00607 neighbors_3d, nb_neighbors, ierror)
00608
00609 if (ierror /= 0) return
00610
00611 #ifdef TODO
00612 case (PSMILe_conserv3D)
00613 call psmile_trs_give_
00614 if (ierror /= 0) return
00615
00616 case (PSMILe_user3D)
00617 call psmile_trs_give_
00618 if (ierror /= 0) return
00619 #endif
00620
00621
00622
00623 case default
00624 ierrp (1) = interpolation_methods (1)
00625 ierror = PRISM_Error_Internal
00626
00627 call psmile_error ( ierror, 'unsupported 3d interpolation method', &
00628 ierrp, 1, __FILE__, __LINE__ )
00629 return
00630 end select
00631
00632
00633
00634
00635
00636
00637
00638
00639 if (global_search) then
00640 if (interpolation_methods(1) /= PSMILe_nnghbr3D) then
00641
00642
00643
00644 call psmile_neigh_global_points (search, extra_search, &
00645 mask_array, mask_shape, src_mask_available, use_how(1), &
00646 grid_id, grid_valid_shape, &
00647 neighbors_3d, ncpl, nb_neighbors, len_cpl, ierror)
00648 if (ierror > 0) return
00649
00650 if (extra_search%n_req > 0) then
00651 call psmile_global_search_dble (comp_info, &
00652 control, len_cpl, &
00653 var_id, grid_valid_shape, search, tgt_grid, &
00654 neighbors_3d, ncpl, nb_neighbors, extra_search, &
00655 interpolation_methods, interpolation_search, &
00656 n_intmethods, &
00657 send_index, src_mask_available, use_mask, &
00658 use_how(1), PRISM_Irrlonlatvrt, ierror)
00659 if (ierror > 0) return
00660 endif
00661
00662 else
00663
00664
00665
00666 endif
00667 endif
00668
00669
00670
00671
00672 if (interpolation_methods(1) /= PSMILe_nnghbr3D) then
00673 call psmile_neigh_extra_points (search, extra_search, &
00674 mask_array, mask_shape, src_mask_available, use_how, &
00675 grid_valid_shape, &
00676 neighbors_3d, ncpl, nb_neighbors, len_cpl, ierror)
00677 if (ierror > 0) return
00678 endif
00679
00680
00681
00682 if (extra_search%n_extra > 0) then
00683
00684
00685
00686
00687 if (global_search) then
00688 call psmile_neigh_extra_search_dble (search, extra_search, &
00689 nb_extra, ierror)
00690 if (ierror > 0) return
00691 endif
00692
00693
00694
00695 if (.not. transformed) then
00696 call psmile_info_trf_coords_3d_dble ( &
00697 mp%coords_pointer%coords_dble(1)%vector, &
00698 mp%coords_pointer%coords_dble(2)%vector, &
00699 mp%coords_pointer%coords_dble(3)%vector, &
00700 mp%coords_pointer%coords_shape, &
00701 grid_valid_shape, &
00702 sinvec, cosvec, ierror)
00703 if (ierror > 0) return
00704
00705 transformed = .true.
00706 endif
00707
00708
00709
00710 do ipart = 1, search%npart
00711 if (len_cpl (ipart) > 0) then
00712
00713 call psmile_neigh_nearestx_3d_dble (grid_id, &
00714 tgt_grid(1)%vector, tgt_grid(2)%vector, &
00715 tgt_grid(3)%vector, &
00716 mp%coords_pointer%coords_dble(1)%vector, &
00717 mp%coords_pointer%coords_dble(2)%vector, &
00718 mp%coords_pointer%coords_dble(3)%vector, &
00719 mp%coords_pointer%coords_shape, &
00720 mask_array, mask_shape, src_mask_available, &
00721 sinvec%vector, cosvec%vector, grid_valid_shape, &
00722 srcloc, ncpl, nprev, len_cpl (ipart), &
00723 neighbors_3d, nb_extra, &
00724 extra_search, ierror)
00725
00726 if (ierror > 0) return
00727
00728 nprev = nprev + len_cpl (ipart)
00729 end if
00730 end do
00731
00732
00733
00734 call psmile_neigh_extra_search_clean (search, extra_search, ierror)
00735 endif
00736
00737
00738
00739
00740 if (transformed) then
00741 Deallocate (cosvec%vector)
00742 Deallocate (sinvec%vector)
00743 endif
00744
00745
00746
00747
00748
00749
00750
00751
00752 Allocate (neighbors(ncpl, nb_neighbors), STAT = ierror)
00753 if ( ierror > 0 ) then
00754 ierrp (1) = ierror
00755 ierrp (2) = ncpl * nb_neighbors
00756
00757 ierror = PRISM_Error_Alloc
00758 call psmile_error ( ierror, 'neighbors', &
00759 ierrp, 2, __FILE__, __LINE__ )
00760 return
00761 endif
00762
00763 call psmile_compact_neighbors_3d (neighbors_3d, ncpl, nb_neighbors, &
00764 grid_valid_shape, extra_search, &
00765 send_info, neighbors, ierror)
00766 if (ierror /= 0) return
00767
00768
00769
00770 if (interpolation_methods(1) /= PSMILe_nnghbr3D) then
00771 call psmile_move0_neighbors (neighbors, ncpl, nb_neighbors, &
00772 ierror)
00773 if (ierror > 0) return
00774 endif
00775
00776 #ifdef DEBUG_TRACE
00777 if (ictl <= ncpl) then
00778 print *, 'i, dstijk (:, i)', ictl, dstijk (:, ictl)
00779
00780 do i = 1, nb_neighbors
00781 print *, 'neighbors_3d (:, ictl)', neighbors_3d (:, ictl, i)
00782 end do
00783
00784 print *, 'neighbors (ictl)', neighbors (ictl, 1:nb_neighbors)
00785 endif
00786 #endif
00787
00788 Deallocate (neighbors_3d)
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805
00806 src_size = send_info%n_list
00807 local_src_size = src_size - send_info%num2recv
00808
00809
00810
00811 if ( src_size == 0 ) then
00812 send_index = PRISM_UNDEFINED
00813
00814 go to 1000
00815 endif
00816
00817 allocate_src = .true.
00818
00819
00820
00821 if (src_mask_available .and. .not. use_mask) then
00822 id_src_mask = 1
00823
00824 Allocate (src_mask(src_size), stat = ierror)
00825 if ( ierror > 0 ) then
00826 ierrp (1) = ierror
00827 ierrp (2) = src_size
00828
00829 ierror = PRISM_Error_Alloc
00830 call psmile_error ( ierror, 'src_mask', &
00831 ierrp, 2, __FILE__, __LINE__ )
00832 return
00833 endif
00834
00835
00836
00837
00838 call psmile_ext_compact_list_log2int (send_info, &
00839 extra_search, mask_array, mask_shape, &
00840 grid_valid_shape, src_mask, src_size, ierror)
00841 if (ierror /= 0) return
00842
00843 else
00844
00845 id_src_mask = 0
00846 src_mask => dummy_mask
00847
00848 endif
00849
00850
00851
00852 if (send_info%send_entire_valid_shape .and. &
00853 send_info%num2recv == 0) then
00854 if (maxval (abs(mp%coords_pointer%coords_shape-grid_valid_shape)) == &
00855 0) then
00856 allocate_src = .false.
00857
00858 nb_corners = 0
00859
00860 call psmile_trs_set_src_epio3d_dble ( &
00861 epio_id, trs_rank, src_size, src_size, &
00862 mp%coords_pointer%coords_dble(1)%vector, &
00863 mp%coords_pointer%coords_dble(2)%vector, &
00864 mp%coords_pointer%coords_dble(3)%vector, &
00865 id_src_mask, src_mask, ierror)
00866 if (ierror /= 0) return
00867
00868 #ifdef DEBUG_TRACE
00869 if (ictl <= ncpl) then
00870
00871 ipart = 1
00872 call psmile_print_3d_coord_dble (coords(1,ipart)%vector, &
00873 coords(2,ipart)%vector, &
00874 coords(3,ipart)%vector, &
00875 shape (:, :, ipart), &
00876 dstijk (:, ictl), 1, 'searched')
00877
00878 print *, 'coords (ictl)'
00879 do i = 1, nb_neighbors
00880 if (neighbors (ictl, i) > 0) then
00881 print *, 'i', i, &
00882 mp%coords_pointer%coords_dble(1)%vector(neighbors (ictl, i)),&
00883 mp%coords_pointer%coords_dble(2)%vector(neighbors (ictl, i)),&
00884 mp%coords_pointer%coords_dble(3)%vector(neighbors (ictl, i))
00885 endif
00886 end do
00887 endif
00888 #endif
00889 endif
00890 endif
00891
00892
00893
00894
00895
00896 if (allocate_src) then
00897 do i = 1, ndim_3d
00898 Allocate (src_grid(i)%vector(src_size), stat = ierror)
00899 if ( ierror > 0 ) then
00900 ierrp (1) = ierror
00901 ierrp (2) = src_size
00902
00903 ierror = PRISM_Error_Alloc
00904 call psmile_error ( ierror, 'src_grid(i)%vector', &
00905 ierrp, 2, __FILE__, __LINE__ )
00906 return
00907 endif
00908
00909 call psmile_ext_compact_list_3d_dble ( &
00910 send_info, &
00911 mp%coords_pointer%coords_dble(i)%vector, &
00912 mp%coords_pointer%coords_shape, &
00913 grid_valid_shape, &
00914 src_grid(i)%vector, src_size, &
00915 ierror)
00916
00917 if (ierror /= 0) return
00918 end do
00919
00920
00921
00922 if (send_info%num2recv > 0) then
00923 do i = 1, ndim_3d
00924 ip = local_src_size
00925
00926 do n = 1, send_info%nrecv
00927 src_grid(i)%vector(ip+1:ip+send_info%len_sent(n)) = &
00928 extra_search%dble_bufs(n)%vector( &
00929 (i-1)*send_info%len_sent(n)+1: &
00930 i *send_info%len_sent(n))
00931
00932 ip = ip + send_info%len_sent (n)
00933 end do
00934 end do
00935
00936 do n = 1, send_info%nrecv
00937 Deallocate (extra_search%dble_bufs(n)%vector)
00938 end do
00939 endif
00940
00941
00942 #ifdef DEBUG_TRACE
00943 if (ictl <= ncpl) then
00944
00945 ipart = 1
00946 call psmile_print_3d_coord_dble (coords(1,ipart)%vector, &
00947 coords(2,ipart)%vector, &
00948 coords(3,ipart)%vector, &
00949 shape (:, :, ipart), &
00950 dstijk (:, ictl), 1, 'searched')
00951
00952 print *, 'coords (ictl)'
00953 do i = 1, nb_neighbors
00954 if (neighbors (ictl, i) > 0) then
00955 print *, 'i', i, &
00956 src_grid(1)%vector(neighbors (ictl, i)),&
00957 src_grid(2)%vector(neighbors (ictl, i)),&
00958 src_grid(3)%vector(neighbors (ictl, i))
00959 endif
00960 end do
00961
00962 endif
00963 #endif
00964
00965 interpolation(1) = interpolation_methods(1)
00966 interpolation(2:ndim_3d) = PSMILe_Undef
00967
00968 call psmile_get_epio_handle (field%comp_id, grid_id, method_id, mask_id, &
00969 interpolation, search%msg_intersections, id_trans_out, &
00970 id_trans_in, search%sender, cpl_id, epio_id, trs_rank, ierror)
00971
00972 if ( epio_id /= PSMILe_undef ) then
00973
00974
00975
00976
00977
00978
00979 call psmile_trs_set_triple_links (id_trans_out, id_trans_in, &
00980 epio_id, trs_rank, ierror)
00981 if (ierror /= 0) return
00982
00983
00984
00985 do i = 1, ndim_3d
00986 Deallocate (src_grid(i)%vector)
00987 end do
00988
00989 if (id_src_mask > 0) then
00990 Deallocate (src_mask)
00991 endif
00992
00993 if (use_target_grid) then
00994 do i = 1, ndim_3d
00995 Deallocate (tgt_grid(i)%vector)
00996 end do
00997 endif
00998
00999 else
01000
01001
01002
01003
01004
01005
01006
01007 nb_corners = 0
01008
01009 call psmile_trs_set_src_epio3d_dble (epio_id, trs_rank, src_size, src_size, &
01010 src_grid(1)%vector, src_grid(2)%vector, src_grid(3)%vector, &
01011 id_src_mask, src_mask, ierror)
01012 if (ierror /= 0) return
01013
01014
01015
01016 do i = 1, ndim_3d
01017 Deallocate (src_grid(i)%vector)
01018 end do
01019 endif
01020
01021 if (id_src_mask == 1) then
01022 Deallocate (src_mask)
01023 endif
01024
01025
01026
01027
01028
01029
01030
01031
01032 id_tgt_mask = 0
01033 nb_corners = 1
01034
01035 call psmile_trs_set_tgt_epio3d_dble ( &
01036 epio_id, trs_rank, ncpl, nb_corners, &
01037 tgt_grid(1)%vector, tgt_grid(2)%vector, tgt_grid(3)%vector, &
01038 id_tgt_mask, dummy_mask, ierror)
01039 if (ierror /= 0) return
01040
01041
01042
01043 if (.not. use_target_grid) then
01044 do i = 1, ndim_3d
01045 Deallocate (tgt_grid(i)%vector)
01046 end do
01047 endif
01048
01049
01050
01051
01052
01053 call psmile_trs_set_triple_links (id_trans_out, id_trans_in, &
01054 epio_id, trs_rank, ierror)
01055 if (ierror /= 0) return
01056
01057
01058
01059 call psmile_trs_give_neighbors3d (epio_id, trs_rank, &
01060 ncpl, nb_neighbors, neighbors, &
01061 ierror)
01062 if (ierror /= 0) return
01063
01064
01065
01066 cpl_list(cpl_id)%epio_id = epio_id
01067 cpl_list(cpl_id)%trs_rank = trs_rank
01068
01069 endif
01070
01071 send_info%epio_id = epio_id
01072 send_info%trs_rank = trs_rank
01073
01074 #ifdef DEBUG
01075 print *, trim(ch_id), 'index, epio_id, n_list ', &
01076 send_index, epio_id, send_info%n_list, &
01077 send_info%send_entire_valid_shape
01078 #endif
01079
01080
01081
01082 1000 continue
01083
01084 call psmile_locations_dealloc (send_info, ierror)
01085
01086
01087
01088 #ifdef VERBOSE
01089 print 9980, trim(ch_id), ierror
01090
01091 call psmile_flushstd
01092 #endif /* VERBOSE */
01093
01094
01095
01096
01097 #ifdef VERBOSE
01098
01099 9990 format (1x, a, ': psmile_info_trs_locs_3d_dble: method_id', i3, &
01100 ' to ', i3, '(', i2, ')')
01101 9980 format (1x, a, ': psmile_info_trs_locs_3d_dble: eof ierror =', i3)
01102
01103 9970 format (1x, a, ': psmile_info_trs_locs_3d_dble: eof src_size=0, ierror =', i3)
01104
01105 #endif /* VERBOSE */
01106
01107 #ifdef DEBUG
01108 9960 format (1x, 'Taskout%In_channel(', i3, '): origin_type', i5, &
01109 ', remote comp id', i4, &
01110 ', global_transi_id', i4, ', remote_transi_id', i4)
01111 #endif
01112
01113 end subroutine PSMILe_info_trs_locs_3d_dble