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