3 subroutine qneg3 (subnam ,idx ,ncol ,ncold ,lver ,lconst_beg , &
5 !-----------------------------------------------------------------------
8 ! Check moisture and tracers for minimum value, reset any below
9 ! minimum value to minimum value and return information to allow
10 ! warning message to be printed. The global average is NOT preserved.
13 ! <Describe the algorithm(s) used in the routine.>
14 ! <Also include any applicable external references.>
18 !-----------------------------------------------------------------------
19 use shr_kind_mod, only: r8 => shr_kind_r8
21 use cam_logfile, only: iulog
23 use module_cam_support, only: iulog
27 !------------------------------Arguments--------------------------------
31 character*(*), intent(in) :: subnam ! name of calling routine
33 integer, intent(in) :: idx ! chunk/latitude index
34 integer, intent(in) :: ncol ! number of atmospheric columns
35 integer, intent(in) :: ncold ! declared number of atmospheric columns
36 integer, intent(in) :: lver ! number of vertical levels in column
37 integer, intent(in) :: lconst_beg ! beginning constituent
38 integer, intent(in) :: lconst_end ! ending constituent
40 real(r8), intent(in) :: qmin(lconst_beg:lconst_end) ! Global minimum constituent concentration
43 ! Input/Output arguments
45 real(r8), intent(inout) :: q(ncold,lver,lconst_beg:lconst_end) ! moisture/tracer field
47 !---------------------------Local workspace-----------------------------
49 integer indx(ncol,lver) ! array of indices of points < qmin
50 integer nval(lver) ! number of points < qmin for 1 level
51 integer nvals ! number of values found < qmin
54 integer i,ii,k ! longitude, level indices
55 integer m ! constituent index
56 integer iw,kw ! i,k indices of worst violator
58 logical found ! true => at least 1 minimum violator found
60 real(r8) worst ! biggest violator
62 !-----------------------------------------------------------------------
65 do m=lconst_beg,lconst_end
71 ! Test all field values for being less than minimum value. Set q = qmin
72 ! for all such points. Trace offenders and identify worst one.
80 if (q(i,k,m) < qmin(m)) then
91 nvals = nvals + nval(k)
93 !cdir nodep,altcode=loopcnt
96 if (q(i,k,m) < worst) then
101 if (iwtmp /= -1 ) kw = k
102 if (iwtmp /= -1 ) iw = indx(iwtmp,k)
103 !cdir nodep,altcode=loopcnt
110 if (found .and. abs(worst)>1.e-12_r8) then
111 write(iulog,9000)subnam,m,idx,nvals,qmin(m),worst,iw,kw
113 call wrf_debug(400,iulog)
119 9000 format(' QNEG3 from ',a,':m=',i3,' lat/lchnk=',i5, &
120 ' Min. mixing ratio violated at ',i4,' points. Reset to ', &
121 1p,e8.1,' Worst =',e8.1,' at i,k=',i4,i3)
126 #if ( defined MODAL_AERO )
127 subroutine qneg3_modalx1 (subnam ,idx ,ncol ,ncold ,lver ,lconst_beg , &
128 lconst_end ,qmin ,q ,qneg3_worst_thresh )
129 !-----------------------------------------------------------------------
132 ! Check moisture and tracers for minimum value, reset any below
133 ! minimum value to minimum value and return information to allow
134 ! warning message to be printed. The global average is NOT preserved.
137 ! <Describe the algorithm(s) used in the routine.>
138 ! <Also include any applicable external references.>
140 ! Author: J. Rosinski
142 !-----------------------------------------------------------------------
143 use shr_kind_mod, only: r8 => shr_kind_r8
145 use cam_logfile, only: iulog
147 use module_cam_support, only: iulog
151 !------------------------------Arguments--------------------------------
155 character*(*), intent(in) :: subnam ! name of calling routine
157 integer, intent(in) :: idx ! chunk/latitude index
158 integer, intent(in) :: ncol ! number of atmospheric columns
159 integer, intent(in) :: ncold ! declared number of atmospheric columns
160 integer, intent(in) :: lver ! number of vertical levels in column
161 integer, intent(in) :: lconst_beg ! beginning constituent
162 integer, intent(in) :: lconst_end ! ending constituent
164 real(r8), intent(in) :: qmin(lconst_beg:lconst_end) ! Global minimum constituent concentration
165 real(r8), intent(in) :: qneg3_worst_thresh(lconst_beg:lconst_end)
166 ! thresholds for reporting violators
169 ! Input/Output arguments
171 real(r8), intent(inout) :: q(ncold,lver,lconst_beg:lconst_end) ! moisture/tracer field
173 !---------------------------Local workspace-----------------------------
175 integer indx(ncol,lver) ! array of indices of points < qmin
176 integer nval(lver) ! number of points < qmin for 1 level
177 integer nvals ! number of values found < qmin
180 integer i,ii,k ! longitude, level indices
181 integer m ! constituent index
182 integer iw,kw ! i,k indices of worst violator
184 logical found ! true => at least 1 minimum violator found
186 real(r8) worst ! biggest violator
187 real(r8) tmp_worst_thresh
189 !-----------------------------------------------------------------------
192 do m=lconst_beg,lconst_end
198 ! Test all field values for being less than minimum value. Set q = qmin
199 ! for all such points. Trace offenders and identify worst one.
207 if (q(i,k,m) < qmin(m)) then
216 if (nval(k) > 0) then
218 nvals = nvals + nval(k)
220 !cdir nodep,altcode=loopcnt
223 if (q(i,k,m) < worst) then
228 if (iwtmp /= -1 ) kw = k
229 if (iwtmp /= -1 ) iw = indx(iwtmp,k)
230 !cdir nodep,altcode=loopcnt
238 tmp_worst_thresh = 1.0e-12_r8
239 if (qneg3_worst_thresh(m) > 0.0_r8) &
240 tmp_worst_thresh = qneg3_worst_thresh(m)
242 if (found .and. abs(worst)>tmp_worst_thresh) then
243 write(iulog,9000)subnam,m,idx,nvals,qmin(m),worst,iw,kw
245 call wrf_debug(400,iulog)
251 9000 format(' QNEG3 from ',a,':m=',i3,' lat/lchnk=',i5, &
252 ' Min. mixing ratio violated at ',i4,' points. Reset to ', &
253 1p,e8.1,' Worst =',e8.1,' at i,k=',i4,i3)
254 end subroutine qneg3_modalx1