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