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