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