00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 subroutine psmile_neigh_tricu_gauss2_irreg ( &
00013 grid_id, grid_valid_shape, interpolation_mode, &
00014 srcloc, srclocz, virtual_cell, nlocs, nloc, nprev, &
00015 virtual_cell_available, neighbors_3d, num_neigh, &
00016 neigh_bascule, ierror)
00017
00018
00019
00020 use PRISM_constants
00021 use PSMILe, dummy_interface => PSMILe_Neigh_tricu_gauss2_irreg
00022 use psmile_grid_reduced_gauss
00023 #ifdef DEBUG_TRACE
00024 use psmile_debug_trace
00025 #endif
00026
00027 implicit none
00028
00029
00030
00031 Integer, Intent (In) :: nloc
00032
00033
00034
00035 Integer, Intent (In) :: nlocs (ndim_2d)
00036
00037
00038
00039 Integer, Intent (In) :: srcloc (ndim_2d, nlocs (1))
00040 Integer, Intent (In) :: srclocz (nlocs (2))
00041
00042
00043
00044 Integer, Intent (In) :: nprev
00045
00046
00047
00048 Integer, Intent (In) :: virtual_cell (nlocs(1))
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063 Integer, Intent (In) :: grid_id
00064
00065
00066
00067 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00068
00069
00070
00071 Integer, Intent (In) :: interpolation_mode
00072
00073
00074
00075
00076
00077
00078
00079 Integer, Intent (In) :: num_neigh
00080
00081
00082
00083 Logical, Intent (In) :: virtual_cell_available
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101 Integer, Intent (Out) :: neighbors_3d (ndim_3d, nloc, num_neigh)
00102
00103
00104
00105 Integer, Intent (Out) :: neigh_bascule(nloc)
00106
00107
00108
00109 Integer, Intent (Out) :: ierror
00110
00111
00112
00113
00114
00115
00116
00117 integer, parameter :: shift_3d(ndim_3d, 9) =
00118 reshape ((/-1,-1,0, 0,-1,0, +1,-1,0,
00119 -1, 0,0, 0, 0,0, +1, 0,0,
00120 -1,+1,0, 0,+1,0, +1,+1,0/),
00121 (/ndim_3d,9/))
00122
00123
00124
00125
00126
00127 Integer :: i, n
00128
00129
00130
00131 Integer, Parameter :: nerrp = 2
00132 Integer :: ierrp (nerrp)
00133
00134 #ifdef DEBUG_TRACE
00135
00136
00137
00138 integer :: ictl_neigh_global_1d(16), ictl_neigh_local_1d(16)
00139 #endif
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165 Character(len=len_cvs_string), save :: mycvs =
00166 '$Id: psmile_neigh_tricu_gauss2_irreg.F90 2968 2011-02-18 15:26:34Z hanke $'
00167
00168
00169
00170 #ifdef VERBOSE
00171 print 9990, trim(ch_id)
00172 call psmile_flushstd
00173 #endif /* VERBOSE */
00174
00175
00176
00177
00178 #ifdef PRISM_ASSERTION
00179
00180
00181
00182 if (interpolation_mode /= 2) then
00183 ierrp (1) = interpolation_mode
00184 ierror = PRISM_Error_Internal
00185
00186 call psmile_error ( ierror, 'unsupported interpolation mode in psmile_neigh_tricu_gauss2_irreg', &
00187 ierrp, 1, __FILE__, __LINE__ )
00188 return
00189 endif
00190
00191
00192
00193 if (num_neigh < 16) then
00194 print 9970, trim(ch_id), num_neigh, 16
00195 call psmile_assert (__FILE__, __LINE__, &
00196 'Number of neighbors is insufficient')
00197 endif
00198
00199
00200
00201 if (nloc < nprev + nlocs (1) * nlocs (2) ) then
00202 print *, trim(ch_id), " nloc, nprev, nlocs ", nloc, nprev, nlocs
00203 call psmile_assert (__FILE__, __LINE__, &
00204 'nloc < nprev + PRODUCT (nlocs) ')
00205 endif
00206
00207
00208
00209
00210 do i = 1, nlocs (1)
00211 if ( abs(srcloc(1,i)) < grid_valid_shape(1,1) .or. &
00212 abs(srcloc(1,i)) > grid_valid_shape(2,1) ) exit
00213 end do
00214
00215 if (i <= nlocs(1)) then
00216 print *, "grid valid shape is ", grid_valid_shape
00217 print *, "Incorrect index in srcloc, i", i, srcloc(1,i)
00218 call psmile_assert (__FILE__, __LINE__, 'wrong index')
00219 endif
00220
00221
00222 do i = 1, nlocs (2)
00223 if ( srclocz(i) < grid_valid_shape(1,3)-1 .or. &
00224 srclocz(i) > grid_valid_shape(2,3)) exit
00225 end do
00226
00227 if (i <= nlocs(2)) then
00228 print *, "Incorrect index in srclocs(2), i", i, srclocz(i), &
00229 grid_valid_shape(1,3)-1, grid_valid_shape(2,3)
00230 call psmile_assert (__FILE__, __LINE__, 'wrong index')
00231 endif
00232 #endif
00233
00234
00235
00236 ierror = 0
00237
00238
00239
00240
00241
00242 if (virtual_cell_available) then
00243
00244 do i = 1, nlocs(1)
00245
00246 n = i + nprev
00247
00248 if (virtual_cell(i) == 0) then
00249
00250 neighbors_3d(:,n,1:16) = psmile_gauss_get_bicu_stencil (grid_id, &
00251 psmile_gauss_local_1d_to_3d (grid_id, abs(srcloc(1,i))))
00252
00253 else
00254
00255 neighbors_3d(:,n,1:16) = psmile_gauss_shift_bicu_stencil (grid_id, &
00256 psmile_gauss_local_1d_to_3d (grid_id, abs(srcloc(1,i))), &
00257 shift_3d(:,virtual_cell(i)))
00258
00259 endif
00260 enddo
00261 else
00262 do i = 1, nlocs(1)
00263
00264 n = i + nprev
00265
00266 neighbors_3d(:,n,1:16) = psmile_gauss_get_bicu_stencil (grid_id, &
00267 psmile_gauss_local_1d_to_3d (grid_id, abs(srcloc(1,i))))
00268 enddo
00269 endif
00270
00271 #ifdef DEBUG_TRACE
00272 if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00273 print *, "neighbors_3d (3d):"
00274 do i = 1, 16
00275 print *, ' neighbors_3d (:,ictl,i) ', i, ':', neighbors_3d (:, ictl, i)
00276 enddo
00277 endif
00278 #endif
00279
00280
00281
00282 where (neighbors_3d(2,nprev+1:nprev+nlocs(1),6) >= grids(grid_id)%nbr_latitudes-1 .or. &
00283 neighbors_3d(2,nprev+1:nprev+nlocs(1),6) == 1)
00284
00285 neigh_bascule(nprev+1:nprev+nlocs(1)) = 1
00286
00287 else where
00288
00289 neigh_bascule(nprev+1:nprev+nlocs(1)) = 0
00290
00291 end where
00292
00293 #ifdef DEBUG_TRACE
00294 if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00295 print *, "neigh_bascule(ictl)", neigh_bascule(ictl)
00296 endif
00297 #endif
00298
00299
00300
00301 do i = 1, 16
00302
00303 neighbors_3d(1,nprev+1:nprev+nlocs(1),i) = &
00304 psmile_gauss_3d_to_global_1d (grid_id, neighbors_3d(:,nprev+1:nprev+nlocs(1),i), nlocs(1))
00305 neighbors_3d(2,nprev+1:nprev+nlocs(1),i) = 1
00306 #ifdef DEBUG_TRACE
00307 if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00308 ictl_neigh_global_1d(i) = neighbors_3d (1, ictl, i)
00309 endif
00310 #endif
00311
00312 neighbors_3d(1,nprev+1:nprev+nlocs(1),i) = &
00313 psmile_gauss_1d_global_to_local (grid_id, neighbors_3d(1,nprev+1:nprev+nlocs(1),i), nlocs(1), 0)
00314
00315 #ifdef DEBUG_TRACE
00316 if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00317 ictl_neigh_local_1d(i) = neighbors_3d (1, ictl, i)
00318 endif
00319 #endif
00320 enddo
00321
00322 #ifdef DEBUG_TRACE
00323 if (ictl > nprev .and. ictl <= nprev + nlocs(1)) then
00324 print *, "neighbors_3d (1d global):"
00325 print *, ' neighbors_3d (1,ictl,1:16) ', ictl_neigh_global_1d
00326 print *, "neighbors_3d (1d local):"
00327 print *, ' neighbors_3d (1,ictl,1:16) ', ictl_neigh_local_1d
00328 endif
00329 #endif
00330
00331 if ( num_neigh == 32 ) then
00332
00333 do n = 1 + nprev, nprev + nlocs(1)
00334 neighbors_3d (1:2, n, 17:32) = neighbors_3d (1:2, n, 1:16)
00335 enddo
00336 endif
00337
00338
00339
00340 do i = 2, nlocs (2)
00341 neighbors_3d (1:2, nprev+(i-1)*nlocs(1)+1:nprev+i*nlocs(1),:) = &
00342 neighbors_3d (1:2, nprev+ 1:nprev+ nlocs(1),:)
00343 end do
00344
00345
00346
00347 do i = 1, nlocs (2)
00348 neighbors_3d (3, nprev+(i-1)*nlocs(1)+1:nprev+i*nlocs(1), :) = srclocz (i)
00349 end do
00350
00351
00352
00353 if (grid_valid_shape(2,3) - grid_valid_shape(1,3) > 0 &
00354 .and. num_neigh == 32) then
00355 do n = 17, 32
00356 neighbors_3d (3, nprev+1:nprev+nlocs(1)*nlocs(2),n) = &
00357 neighbors_3d (3, nprev+1:nprev+nlocs(1)*nlocs(2), 1) + 1
00358 enddo
00359 endif
00360
00361
00362
00363 #ifdef VERBOSE
00364 print 9980, trim(ch_id)
00365
00366 call psmile_flushstd
00367 #endif /* VERBOSE */
00368
00369
00370
00371 9990 format (1x, a, ': psmile_neigh_tricu_gauss2_irreg')
00372 9980 format (1x, a, ': psmile_neigh_tricu_gauss2_irreg: eof')
00373
00374 #ifdef PRISM_ASSERTION
00375 9970 format (1x, a, ': number of neighbors', i2, '; expected at least', i2)
00376 #endif
00377
00378 end subroutine psmile_neigh_tricu_gauss2_irreg