00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_bbcells_1d_dble (array, shape, range, &
00012 corners, corner_shape, nbr_corners, &
00013 chmin, chmax, levdim, &
00014 ierror)
00015
00016
00017
00018 use PRISM_constants
00019
00020 use PSMILe, dummy_interface => PSMILe_Bbcells_1d_dble
00021
00022 implicit none
00023
00024
00025
00026 Integer, Intent (In) :: shape (2)
00027
00028
00029
00030 Double Precision, Intent (In) :: array(shape(1):shape(2))
00031
00032
00033
00034 Integer, Intent (In) :: range (2)
00035
00036
00037
00038
00039 Integer, Intent (In) :: corner_shape (2)
00040
00041
00042
00043
00044 Integer, Intent (In) :: nbr_corners
00045
00046
00047
00048
00049 Double Precision, Intent (In) :: corners (corner_shape(1):
00050 corner_shape(2),
00051 nbr_corners)
00052
00053
00054
00055 Integer, Intent (In) :: levdim
00056
00057
00058
00059
00060
00061 Double Precision, Intent (Out) :: chmin (range(1)-1:range(2))
00062
00063
00064
00065 Double Precision, Intent (Out) :: chmax (range(1)-1:range(2))
00066
00067
00068
00069 Integer, Intent (Out) :: ierror
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079 Double Precision, Parameter :: epsp1 = 1.0d0 + 5d-6
00080
00081
00082
00083 Integer :: i, ibeg, iend
00084
00085 Double Precision :: chmin_corner (2),
00086 chmax_corner (2)
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118 Character(len=len_cvs_string), save :: mycvs =
00119 '$Id: psmile_bbcells_1d_dble.F90 2807 2010-12-08 09:56:19Z redler $'
00120
00121
00122
00123 #ifdef VERBOSE
00124 print 9990, trim(ch_id)
00125
00126 call psmile_flushstd
00127 #endif /* VERBOSE */
00128
00129 #ifdef PRISM_ASSERTION
00130 if (range(1) < corner_shape(1) .or. range(2) > corner_shape(2)) then
00131 print *, trim(ch_id), ": PSMILe_bbcells_1d_dble: range, corner", &
00132 range, corner_shape
00133 call psmile_assert (__FILE__, __LINE__, &
00134 "levdim_corner must be >= range(2)-range(1)")
00135 endif
00136
00137 if (levdim < range(2)-range(1)+2) then
00138 print *, trim(ch_id), ": PSMILe_bbcells_1d_dble: levdim, range", levdim
00139 call psmile_assert (__FILE__, __LINE__, &
00140 "levdim < range(2)-range(1)+2")
00141 endif
00142
00143 if (nbr_corners /= 2) then
00144 print *, trim(ch_id), ": PSMILe_bbcells_1d_dble: nbr_corner", &
00145 nbr_corners
00146 call psmile_assert (__FILE__, __LINE__, &
00147 "nbr_corners != 2")
00148 endif
00149 #endif
00150
00151
00152
00153 ierror = 0
00154
00155 ibeg = range (1)
00156 iend = range (2)
00157
00158
00159
00160 chmin_corner (1) = min (corners(ibeg, 1), corners(ibeg, 2))
00161 chmax_corner (1) = max (corners(ibeg, 1), corners(ibeg, 2))
00162
00163 chmin_corner (2) = min (corners(iend, 1), corners(iend, 2))
00164 chmax_corner (2) = max (corners(iend, 1), corners(iend, 2))
00165
00166 if (ibeg == iend) then
00167
00168 #ifdef PRISM_ASSERTION
00169 if (chmin_corner (1) > array (ibeg) .or. &
00170 chmax_corner (1) < array (ibeg)) then
00171 print *, trim(ch_id), ": PSMILe_bbcells_1d_dble: corner", &
00172 chmin_corner (1), chmax_corner (1), array (ibeg)
00173 call psmile_assert (__FILE__, __LINE__, &
00174 "method point should be in corner control volume")
00175 endif
00176 #endif
00177
00178 chmin (ibeg-1) = chmin_corner (1)
00179 chmax (ibeg-1) = array (ibeg)
00180
00181 chmin (iend) = array (iend)
00182 chmax (iend) = chmax_corner (1)
00183
00184 else if ( ibeg < iend) then
00185
00186
00187
00188
00189 do i = ibeg, iend-1
00190 chmin (i) = min (array(i), array (i+1))
00191 enddo
00192
00193
00194 do i = ibeg, iend-1
00195 chmax (i) = max (array(i), array (i+1))
00196 enddo
00197
00198
00199
00200
00201
00202
00203 if (array(ibeg) == chmin (ibeg)) then
00204 chmin (ibeg-1) = chmin_corner (1)
00205 chmax (ibeg-1) = array (ibeg)
00206 else
00207 chmin (ibeg-1) = array (ibeg)
00208 chmax (ibeg-1) = chmax_corner (1)
00209 endif
00210
00211 if (array (iend) == chmax (iend-1) ) then
00212 chmin (iend) = array (iend)
00213 chmax (iend) = chmax_corner (2)
00214 else
00215 chmin (iend) = chmin_corner (2)
00216 chmax (iend) = array (iend)
00217 endif
00218
00219 endif
00220
00221
00222
00223 #ifdef VERBOSE
00224 print 9980, trim(ch_id), ierror
00225
00226 call psmile_flushstd
00227 #endif /* VERBOSE */
00228
00229
00230
00231 9990 format (1x, a, ': psmile_bbcells_1d_dble')
00232 9980 format (1x, a, ': psmile_bbcells_1d_dble: eof ierror =', i3)
00233
00234 end subroutine PSMILe_Bbcells_1d_dble