00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_remove_intersect (inter, idl, idg, npart, &
00012 local_extent_info, global_extent_info, ierror)
00013
00014
00015
00016 use PRISM_constants
00017
00018 use PSMILe, dummy_interface => PSMILe_Remove_intersect
00019
00020 implicit none
00021
00022
00023
00024 Type (Enddef_extent_info), Intent (In) ::
00025 local_extent_info (*)
00026
00027
00028
00029
00030 Type (Enddef_extent_info), Intent (In) ::
00031 global_extent_info (*)
00032
00033
00034
00035
00036
00037
00038 Integer, Intent (InOut) :: npart
00039
00040
00041
00042 Real (PSMILe_float_kind), Intent (InOut) :: inter (2, ndim_3d, npart)
00043
00044
00045
00046 Integer, Intent (InOut) :: idl (npart)
00047
00048
00049
00050
00051 Integer, Intent (InOut) :: idg (npart)
00052
00053
00054
00055
00056
00057
00058 Integer, Intent (Out) :: ierror
00059
00060
00061
00062
00063
00064
00065
00066 Double Precision :: gap = 0.0d0
00067
00068
00069
00070 Integer :: i, k, l, lbeg
00071 Integer :: nfit
00072 Logical :: fit (ndim_3d)
00073
00074 Integer :: trans_diff (npart)
00075 Real (PSMILe_float_kind) :: common (2, ndim_3d, npart)
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094 Character(len=len_cvs_string), save :: mycvs =
00095 '$Id: psmile_remove_intersect.F90 3248 2011-06-23 13:03:19Z coquart $'
00096
00097
00098
00099
00100
00101 #ifdef VERBOSE
00102 print 9990, trim(ch_id), npart
00103
00104 call psmile_flushstd
00105
00106 #endif /* VERBOSE */
00107
00108 ierror = 0
00109
00110 do l = 1, npart
00111 trans_diff (l) = local_extent_info(idl(l))%tr_code - &
00112 global_extent_info(idg(l))%tr_code
00113 end do
00114
00115
00116
00117 k = 1
00118 do while (k <= npart - 1)
00119
00120
00121
00122 do l = k+1, npart
00123 common (1, :, l) = max (inter (1,:,l), inter (1,:,k))
00124 common (2, :, l) = min (inter (2,:,l), inter (2,:,k))
00125 end do
00126
00127
00128
00129
00130
00131 do l = k+1, npart
00132 if (trans_diff (l) /= trans_diff (k)) &
00133 common (1, 1, l) = common (2, 1, l) + 1.0
00134 end do
00135
00136 lbeg = k + 1
00137
00138 do while (lbeg <= npart )
00139 do l = lbeg, npart
00140 if (minval (common (2, :, l) - common (1, :, l)) >= 0.0d0) exit
00141 end do
00142
00143 if (l > npart) exit
00144 lbeg = l + 1
00145
00146 if (maxval (common (2, :, l) - common (1, :, l)) == 0.0d0) cycle
00147
00148
00149
00150
00151
00152
00153
00154 do i = 1, ndim_3d
00155 fit (i) = maxval (abs (common (1:2,i,l) - inter (1:2,i,l))) &
00156 == 0.0d0
00157 end do
00158 nfit = count (fit)
00159
00160
00161
00162 if (nfit == ndim_3d) then
00163 if (l /= npart) then
00164 inter (:, :, l) = inter (:, :, npart)
00165 common (:, :, l) = common (:, :, npart)
00166 idl (l) = idl (npart)
00167 idg (l) = idg (npart)
00168 trans_diff (l) = trans_diff (npart)
00169 endif
00170
00171 npart = npart - 1
00172 lbeg = l
00173
00174 cycle
00175 endif
00176
00177
00178
00179 if (maxval (abs (common (:,:,l) - inter (:,:,k))) == 0.0d0) then
00180 inter (:, :, k) = inter (:, :, npart)
00181 trans_diff (k) = trans_diff (npart)
00182
00183 idl (k) = idl (npart)
00184 idg (k) = idg (npart)
00185
00186 npart = npart - 1
00187 k = k - 1
00188
00189 exit
00190 endif
00191
00192
00193
00194
00195
00196
00197 if (nfit >= ndim_2d) then
00198 if (fit(1) .and. fit (2)) then
00199
00200 if (common (1, 3, l) == inter (1, 3, l)) then
00201 inter (1, 3, l) = common (2, 3, l) + gap
00202
00203
00204 else if (common (2, 3, l) == inter (2, 3, l)) then
00205
00206 inter (2, 3, l) = common (1, 3, l) - gap
00207
00208 endif
00209
00210 else if (fit (2) .and. fit (3)) then
00211
00212 if (common (1, 1, l) == inter (1, 1, l)) then
00213 inter (1, 1, l) = common (2, 1, l) + gap
00214
00215
00216 else if (common (2, 1, l) == inter (2, 1, l)) then
00217
00218 inter (2, 1, l) = common (1, 1, l) - gap
00219
00220 endif
00221
00222 else
00223
00224 if (common (1, 2, l) == inter (1, 2, l)) then
00225 inter (1, 2, l) = common (2, 2, l) + gap
00226
00227
00228 else if (common (2, 1, l) == inter (2, 1, l)) then
00229
00230 inter (2, 2, l) = common (1, 2, l) - gap
00231
00232 endif
00233 endif
00234 endif
00235
00236 #ifdef WIRKLICH
00237
00238
00239
00240
00241
00242
00243
00244
00245 if (
00246
00247 do i = 1, ndim_3d
00248 fitk (i) = maxval (abs (common (1:2,i,l) - inter (1:2,i,k))) &
00249 == 0
00250 end do
00251 nfitk = count (fitk)
00252
00253 if (nfitk >= ndim_2d) then
00254 if (fitk(1) .and. fitk (2)) then
00255
00256 if (common (1, 3, l) == inter (1, 3, k)) then
00257 inter (1, 3, k) = common (2, 3, l) + gap
00258 k = k - 1
00259 exit
00260
00261 else if (common (2, 3, l) == inter (2, 3, k)) then
00262
00263 inter (2, 3, k) = common (1, 3, l) - gap
00264 k = k - 1
00265 exit
00266
00267 endif
00268
00269 else if (fitk (2) .and. fitk (3)) then
00270
00271 if (common (1, 1, l) == inter (1, 1, k)) then
00272 inter (1, 1, k) = common (2, 1, l) + gap
00273 k = k - 1
00274 exit
00275
00276 else if (common (2, 1, l) == inter (2, 1, k)) then
00277
00278 inter (2, 1, k) = common (1, 1, l) - gap
00279 k = k - 1
00280 exit
00281 endif
00282
00283 else
00284
00285 if (common (1, 2, l) == inter (1, 2, k)) then
00286 inter (1, 2, k) = common (2, 2, l) + gap
00287 k = k - 1
00288 exit
00289
00290 else if (common (2, 1, l) == inter (2, 1, k)) then
00291
00292 inter (2, 2, k) = common (1, 2, l) - gap
00293 k = k - 1
00294 exit
00295 endif
00296 endif
00297 endif
00298 endif
00299 #endif
00300
00301 end do
00302
00303 k = k + 1
00304 end do
00305
00306
00307
00308 #ifdef VERBOSE
00309 print 9980, trim(ch_id), ierror, npart
00310 call psmile_flushstd
00311 #endif /* VERBOSE */
00312
00313
00314
00315
00316 #ifdef VERBOSE
00317
00318 9990 format (1x, a, ': psmile_remove_intersect: npart =', i4)
00319 9980 format (1x, a, ': psmile_remove_intersect: eof ierror =', i3, &
00320 '; npart =', i4)
00321
00322 #endif /* VERBOSE */
00323
00324 end subroutine PSMILe_Remove_intersect