updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / phys / module_cam_mp_qneg3.F
blob2fd53d6a39dec2223efa6a9332a219c118eda622
1 #define WRF_PORT
2 #define MODAL_AERO
3 subroutine qneg3 (subnam  ,idx     ,ncol    ,ncold   ,lver    ,lconst_beg  , &
4                   lconst_end       ,qmin    ,q       )
5 !----------------------------------------------------------------------- 
6
7 ! Purpose: 
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.
11
12 ! Method: 
13 ! <Describe the algorithm(s) used in the routine.> 
14 ! <Also include any applicable external references.> 
15
16 ! Author: J. Rosinski
17
18 !-----------------------------------------------------------------------
19    use shr_kind_mod, only: r8 => shr_kind_r8
20 #ifndef WRF_PORT
21    use cam_logfile,  only: iulog
22 #else
23    use module_cam_support, only: iulog
24 #endif
25    implicit none
27 !------------------------------Arguments--------------------------------
29 ! Input 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
52    integer nn
53    integer iwtmp
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
66       nvals = 0
67       found = .false.
68       worst = 1.e35_r8
69       iw = -1
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.
74 !DIR$ preferstream
75       do k=1,lver
76          nval(k) = 0
77 !DIR$ prefervector
78          nn = 0
79          do i=1,ncol
80             if (q(i,k,m) < qmin(m)) then
81                nn = nn + 1
82                indx(nn,k) = i
83             end if
84          end do
85          nval(k) = nn
86       end do
88       do k=1,lver
89          if (nval(k) > 0) then
90             found = .true.
91             nvals = nvals + nval(k)
92             iwtmp = -1
93 !cdir nodep,altcode=loopcnt
94             do ii=1,nval(k)
95                i = indx(ii,k)
96                if (q(i,k,m) < worst) then
97                   worst = q(i,k,m)
98                   iwtmp = ii
99                end if
100             end do
101             if (iwtmp /= -1 ) kw = k
102             if (iwtmp /= -1 ) iw = indx(iwtmp,k)
103 !cdir nodep,altcode=loopcnt
104             do ii=1,nval(k)
105                i = indx(ii,k)
106                q(i,k,m) = qmin(m)
107             end do
108          end if
109       end do
110       if (found .and. abs(worst)>1.e-12_r8) then
111          write(iulog,9000)subnam,m,idx,nvals,qmin(m),worst,iw,kw
112 #ifdef WRF_PORT
113          call wrf_debug(400,iulog)
114 #endif
115       end if
116    end do
118    return
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)
122 end subroutine qneg3
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 !----------------------------------------------------------------------- 
131 ! Purpose: 
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.
136 ! Method: 
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
144 #ifndef WRF_PORT
145    use cam_logfile,  only: iulog
146 #else
147    use module_cam_support, only: iulog
148 #endif
149    implicit none
151 !------------------------------Arguments--------------------------------
153 ! Input 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
178    integer nn
179    integer iwtmp
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
193       nvals = 0
194       found = .false.
195       worst = 1.e35_r8
196       iw = -1
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.
201 !DIR$ preferstream
202       do k=1,lver
203          nval(k) = 0
204 !DIR$ prefervector
205          nn = 0
206          do i=1,ncol
207             if (q(i,k,m) < qmin(m)) then
208                nn = nn + 1
209                indx(nn,k) = i
210             end if
211          end do
212          nval(k) = nn
213       end do
215       do k=1,lver
216          if (nval(k) > 0) then
217             found = .true.
218             nvals = nvals + nval(k)
219             iwtmp = -1
220 !cdir nodep,altcode=loopcnt
221             do ii=1,nval(k)
222                i = indx(ii,k)
223                if (q(i,k,m) < worst) then
224                   worst = q(i,k,m)
225                   iwtmp = ii
226                end if
227             end do
228             if (iwtmp /= -1 ) kw = k
229             if (iwtmp /= -1 ) iw = indx(iwtmp,k)
230 !cdir nodep,altcode=loopcnt
231             do ii=1,nval(k)
232                i = indx(ii,k)
233                q(i,k,m) = qmin(m)
234             end do
235          end if
236       end do
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
244 #ifdef WRF_PORT
245          call wrf_debug(400,iulog)
246 #endif
247       end if
248    end do
250    return
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
255 #endif