00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_get_faces_irreg2_real (search, extra_search, &
00012 corners1, corners2, corners3, &
00013 corner_shape, nbr_corners, grid_valid_shape, &
00014 neighbors_3d, nloc, num_neigh, &
00015 faces, n_faces, nreq, ierror)
00016
00017
00018
00019 use PRISM
00020
00021 use PSMILe, dummy_interface => PSMILe_Get_faces_irreg2_real
00022
00023 Implicit none
00024
00025
00026
00027 Type (Enddef_search), Intent (In) :: search
00028
00029
00030
00031 Type (Extra_search_info), Intent (In) :: extra_search
00032
00033
00034
00035
00036 Integer, Intent (In) :: corner_shape (2, ndim_3d)
00037
00038
00039
00040 Integer, Intent (In) :: nbr_corners
00041
00042
00043
00044
00045 Integer, Intent (In) :: grid_valid_shape (2, ndim_3d)
00046
00047
00048
00049 Real, Intent (In) :: corners1 (
00050 corner_shape(1,1):corner_shape(2,1),
00051 corner_shape(1,2):corner_shape(2,2),
00052 nbr_corners/2)
00053
00054 Real, Intent (In) :: corners2 (
00055 corner_shape(1,1):corner_shape(2,1),
00056 corner_shape(1,2):corner_shape(2,2),
00057 nbr_corners/2)
00058
00059 Real, Intent (In) :: corners3 (
00060 corner_shape(1,3):corner_shape(2,3),
00061 2)
00062
00063
00064
00065
00066
00067 Integer, Intent (In) :: nloc
00068
00069
00070
00071
00072
00073
00074 Integer, Intent (In) :: num_neigh
00075
00076
00077
00078
00079 Integer, Intent (In) :: neighbors_3d (ndim_3d, nloc,
00080 num_neigh)
00081
00082
00083
00084 Integer, Intent (In) :: nreq
00085
00086
00087
00088
00089
00090
00091 Real, Intent (Out) :: faces (2, ndim_3d, nreq)
00092
00093
00094 Integer, Intent (Out) :: n_faces (nreq)
00095
00096
00097
00098 Integer, Intent (Out) :: ierror
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108 Integer :: i, j, n
00109
00110
00111
00112 Integer :: ipart, nprev
00113
00114 Integer, Pointer :: indices_req (:)
00115 Integer, Pointer :: required (:)
00116 Integer, Pointer :: len_req (:)
00117
00118
00119
00120 Integer :: code
00121
00122 Integer :: ind (ndim_3d)
00123 Real :: box (2, ndim_3d)
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150 Character(len=len_cvs_string), save :: mycvs =
00151 '$Id: psmile_get_faces_irreg2_real.F90 2082 2009-10-21 13:31:19Z hanke $'
00152
00153
00154
00155
00156
00157 #ifdef VERBOSE
00158 print 9990, trim(ch_id)
00159
00160 call psmile_flushstd
00161 #endif /* VERBOSE */
00162
00163 #ifdef PRISM_ASSERTION
00164 #endif
00165
00166 ierror = 0
00167
00168 len_req => extra_search%len_req
00169
00170
00171 n_faces (:) = 1
00172
00173 faces (1, 1, :) = 1.0
00174 faces (2, 1, :) = 0.0
00175
00176
00177
00178
00179
00180
00181
00182 nprev = 0
00183 do ipart = 1, search%npart
00184 if (len_req (ipart) > 0) then
00185 indices_req => extra_search%indices_req(ipart)%vector
00186 required => extra_search%required(ipart)%vector
00187
00188
00189
00190 code = 1
00191
00192 do n = 1, num_neigh
00193
00194 do j = 1, len_req(ipart)
00195 if (IAND (required (j), code) == 0) cycle
00196
00197 i = indices_req (j)
00198
00199 if (neighbors_3d (1,i,n) < grid_valid_shape(1,1)) then
00200 ind (1) = grid_valid_shape(1,1)
00201 else
00202 ind (1) = min (grid_valid_shape(2,1), neighbors_3d (1,i,n))
00203 endif
00204
00205 if (neighbors_3d (2,i,n) < grid_valid_shape(1,2)) then
00206 ind (2) = grid_valid_shape (1,2)
00207 else
00208 ind (2) = min (grid_valid_shape(2,2), neighbors_3d (2,i,n))
00209 endif
00210
00211 if (neighbors_3d (3,i,n) < grid_valid_shape(1,3)) then
00212 ind (3) = grid_valid_shape (1,3)
00213 else
00214 ind (3) = min (grid_valid_shape(2,3), neighbors_3d (3,i,n))
00215 endif
00216
00217
00218
00219 box (1, 1) = minval(corners1 (ind (1), ind(2), :))
00220 box (2, 1) = maxval(corners1 (ind (1), ind(2), :))
00221
00222 box (1, 2) = minval(corners2 (ind (1), ind(2), :))
00223 box (2, 2) = maxval(corners2 (ind (1), ind(2), :))
00224
00225 box (1, 3) = minval(corners3 (ind (3), :))
00226 box (2, 3) = maxval(corners3 (ind (3), :))
00227
00228 if (faces (1, 1, nprev+j) > faces (2, 1, nprev+j)) then
00229
00230
00231
00232 faces (:, :, nprev+j) = box (:, :)
00233
00234 else
00235
00236 faces (1,1, nprev+j) = min (faces (1,1, nprev+j), box (1,1))
00237 faces (2,1, nprev+j) = max (faces (2,1, nprev+j), box (2,1))
00238 faces (1,2, nprev+j) = min (faces (1,2, nprev+j), box (1,2))
00239 faces (2,2, nprev+j) = max (faces (2,2, nprev+j), box (2,2))
00240 faces (1,3, nprev+j) = min (faces (1,3, nprev+j), box (1,3))
00241 faces (2,3, nprev+j) = max (faces (2,3, nprev+j), box (2,3))
00242
00243 endif
00244
00245 end do
00246
00247 code = code * 2
00248 end do
00249 nprev = nprev + len_req (ipart)
00250 endif
00251 end do
00252
00253 #ifdef PRISM_ASSERTION
00254 do j = 1, nreq
00255 if (faces (1,1,j) > faces (2,1,j)) exit
00256 end do
00257
00258 if (j < nreq) then
00259 nprev = 0
00260 do ipart = 1, search%npart
00261 if (nprev + len_req(ipart) >= j) exit
00262 nprev = nprev + len_req(ipart)
00263 end do
00264
00265 indices_req => extra_search%indices_req(ipart)%vector
00266 i = indices_req (j)
00267 print *, 'j, nreq', j, nreq, grid_valid_shape
00268 do n = 1, num_neigh
00269 print *, 'n, ng', n, &
00270 neighbors_3d (1,j,n), &
00271 neighbors_3d (2,j,n), &
00272 neighbors_3d (3,j,n)
00273 end do
00274
00275 call psmile_assert (__FILE__, __LINE__, &
00276 "neighbour info is not consistent")
00277 endif
00278 #endif
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291 #ifdef VERBOSE
00292 print 9980, trim(ch_id), ierror
00293
00294 call psmile_flushstd
00295 #endif /* VERBOSE */
00296
00297
00298
00299
00300 #ifdef VERBOSE
00301
00302 9990 format (1x, a, ': psmile_get_faces_irreg2_real:')
00303 9980 format (1x, a, ': psmile_get_faces_irreg2_real: eof ierror =', i3)
00304
00305 #endif /* VERBOSE */
00306
00307 #ifdef DEBUG
00308 #endif
00309
00310 end subroutine PSMILe_Get_faces_irreg2_real