Line data Source code
1 : !===============================================================================
2 : ! cloud cover output
3 : !===============================================================================
4 : module cloud_cover_diags
5 :
6 : use shr_kind_mod, only: r8=>shr_kind_r8, shr_kind_CS
7 : use ppgrid, only: pcols, pver,pverp
8 : use cam_history, only: addfld, add_default, outfld, horiz_only
9 : use phys_control, only: phys_getopts
10 :
11 : implicit none
12 :
13 : private
14 :
15 : public :: cloud_cover_diags_init
16 : public :: cloud_cover_diags_out
17 :
18 : real(r8) plowmax ! Max prs for low cloud cover range
19 : real(r8) plowmin ! Min prs for low cloud cover range
20 : real(r8) pmedmax ! Max prs for mid cloud cover range
21 : real(r8) pmedmin ! Min prs for mid cloud cover range
22 : real(r8) phghmax ! Max prs for hgh cloud cover range
23 : real(r8) phghmin ! Min prs for hgh cloud cover range
24 : !
25 : parameter (plowmax = 120000._r8,plowmin = 70000._r8, &
26 : pmedmax = 70000._r8,pmedmin = 40000._r8, &
27 : phghmax = 40000._r8,phghmin = 5000._r8)
28 :
29 :
30 : contains
31 :
32 : !===============================================================================
33 : !===============================================================================
34 3072 : subroutine cloud_cover_diags_init(sampling_seq)
35 :
36 : character(len=*), intent(in) :: sampling_seq
37 : logical :: history_amwg ! output the variables used by the AMWG diag package
38 : character(len=shr_kind_CS) :: long_name_string
39 :
40 3072 : call addfld ('CLOUD', (/ 'lev' /), 'A','fraction','Cloud fraction' , sampling_seq=sampling_seq)
41 1536 : call addfld ('CLDTOT',horiz_only, 'A','fraction','Vertically-integrated total cloud' , sampling_seq=sampling_seq)
42 :
43 1536 : write(long_name_string,999) 'Vertically-integrated low cloud from ', plowmin, ' to ', plowmax, ' Pa'
44 1536 : call addfld ('CLDLOW',horiz_only, 'A','fraction',long_name_string , sampling_seq=sampling_seq)
45 :
46 1536 : write(long_name_string,999) 'Vertically-integrated mid-level cloud from ', pmedmin, ' to ', pmedmax, ' Pa'
47 1536 : call addfld ('CLDMED',horiz_only, 'A','fraction',long_name_string , sampling_seq=sampling_seq)
48 :
49 1536 : write(long_name_string,999) 'Vertically-integrated high cloud from ', phghmin, ' to ', phghmax, ' Pa'
50 1536 : call addfld ('CLDHGH',horiz_only, 'A','fraction',long_name_string , sampling_seq=sampling_seq)
51 :
52 : 999 format(A,F7.0,A,F7.0,A)
53 : ! determine the add_default fields
54 1536 : call phys_getopts(history_amwg_out = history_amwg )
55 :
56 1536 : if (history_amwg) then
57 1536 : call add_default ('CLOUD ', 1, ' ')
58 1536 : call add_default ('CLDTOT ', 1, ' ')
59 1536 : call add_default ('CLDLOW ', 1, ' ')
60 1536 : call add_default ('CLDMED ', 1, ' ')
61 1536 : call add_default ('CLDHGH ', 1, ' ')
62 : endif
63 :
64 :
65 1536 : end subroutine cloud_cover_diags_init
66 :
67 : !===============================================================================
68 : !===============================================================================
69 749232 : subroutine cloud_cover_diags_out(lchnk, ncol, cld, pmid, nmxrgn, pmxrgn )
70 :
71 : integer, intent(in) :: lchnk, ncol
72 : real(r8), intent(in) :: cld(pcols,pver)
73 : real(r8), intent(in) :: pmid(pcols,pver)
74 : integer, intent(in) :: nmxrgn(pcols)
75 : real(r8), intent(in) :: pmxrgn(pcols,pverp)
76 :
77 : real(r8) :: cltot(pcols) ! Diagnostic total cloud cover
78 : real(r8) :: cllow(pcols) ! " low cloud cover
79 : real(r8) :: clmed(pcols) ! " mid cloud cover
80 : real(r8) :: clhgh(pcols) ! " hgh cloud cover
81 :
82 749232 : call cldsav (lchnk, ncol, cld, pmid, cltot, cllow, clmed, clhgh, nmxrgn, pmxrgn)
83 :
84 : !
85 : ! Dump cloud field information to history tape buffer (diagnostics)
86 : !
87 749232 : call outfld('CLDTOT ',cltot ,pcols,lchnk)
88 749232 : call outfld('CLDLOW ',cllow ,pcols,lchnk)
89 749232 : call outfld('CLDMED ',clmed ,pcols,lchnk)
90 749232 : call outfld('CLDHGH ',clhgh ,pcols,lchnk)
91 :
92 749232 : call outfld('CLOUD ',cld ,pcols,lchnk)
93 :
94 749232 : end subroutine cloud_cover_diags_out
95 :
96 : !===============================================================================
97 : !===============================================================================
98 749232 : subroutine cldsav(lchnk ,ncol , &
99 : cld ,pmid ,cldtot ,cldlow ,cldmed , &
100 : cldhgh ,nmxrgn ,pmxrgn )
101 : !-----------------------------------------------------------------------
102 : !
103 : ! Purpose:
104 : ! Compute total & 3 levels of cloud fraction assuming maximum-random overlap.
105 : ! Pressure ranges for the 3 cloud levels are specified.
106 : !
107 : ! Method:
108 : ! <Describe the algorithm(s) used in the routine.>
109 : ! <Also include any applicable external references.>
110 : !
111 : ! Author: W. Collins
112 : !
113 : !-----------------------------------------------------------------------
114 :
115 : implicit none
116 : !
117 : !------------------------------Arguments--------------------------------
118 : !
119 : ! Input arguments
120 : !
121 : integer, intent(in) :: lchnk ! chunk identifier
122 : integer, intent(in) :: ncol ! number of atmospheric columns
123 :
124 : real(r8), intent(in) :: cld(pcols,pver) ! Cloud fraction
125 : real(r8), intent(in) :: pmid(pcols,pver) ! Level pressures
126 : real(r8), intent(in) :: pmxrgn(pcols,pverp) ! Maximum values of pressure for each
127 : ! maximally overlapped region.
128 : ! 0->pmxrgn(i,1) is range of pressure for
129 : ! 1st region,pmxrgn(i,1)->pmxrgn(i,2) for
130 : ! 2nd region, etc
131 :
132 : integer, intent(in) :: nmxrgn(pcols) ! Number of maximally overlapped regions
133 : !
134 : ! Output arguments
135 : !
136 : real(r8), intent(out) :: cldtot(pcols) ! Total random overlap cloud cover
137 : real(r8), intent(out) :: cldlow(pcols) ! Low random overlap cloud cover
138 : real(r8), intent(out) :: cldmed(pcols) ! Middle random overlap cloud cover
139 : real(r8), intent(out) :: cldhgh(pcols) ! High random overlap cloud cover
140 :
141 : !
142 : !---------------------------Local workspace-----------------------------
143 : !
144 : integer i,k ! Longitude,level indices
145 : integer irgn(pcols) ! Max-overlap region index
146 : integer max_nmxrgn ! maximum value of nmxrgn over columns
147 : integer ityp ! Type counter
148 : real(r8) clrsky(pcols) ! Max-random clear sky fraction
149 : real(r8) clrskymax(pcols) ! Maximum overlap clear sky fraction
150 : !------------------------------Parameters-------------------------------
151 : real(r8) ptypmin(4)
152 : real(r8) ptypmax(4)
153 :
154 : data ptypmin /phghmin, plowmin, pmedmin, phghmin/
155 : data ptypmax /plowmax, plowmax, pmedmax, phghmax/
156 : !
157 : !-----------------------------------------------------------------------
158 : !
159 : ! Initialize region number
160 : !
161 749232 : max_nmxrgn = -1
162 12510432 : do i=1,ncol
163 12510432 : max_nmxrgn = max(max_nmxrgn,nmxrgn(i))
164 : end do
165 :
166 3746160 : do ityp = 1, 4
167 53038656 : irgn(1:ncol) = 1
168 6588996 : do k =1,max_nmxrgn-1
169 63029744 : do i=1,ncol
170 60032816 : if (pmxrgn(i,irgn(i)) < ptypmin(ityp) .and. irgn(i) < nmxrgn(i)) then
171 0 : irgn(i) = irgn(i) + 1
172 : end if
173 : end do
174 : end do
175 : !
176 : ! Compute cloud amount by estimating clear-sky amounts
177 : !
178 50041728 : clrsky(1:ncol) = 1.0_r8
179 50041728 : clrskymax(1:ncol) = 1.0_r8
180 80917056 : do k = 1, pver
181 1304081856 : do i=1,ncol
182 1301084928 : if (pmid(i,k) >= ptypmin(ityp) .and. pmid(i,k) <= ptypmax(ityp)) then
183 493970400 : if (pmxrgn(i,irgn(i)) < pmid(i,k) .and. irgn(i) < nmxrgn(i)) then
184 695062 : irgn(i) = irgn(i) + 1
185 695062 : clrsky(i) = clrsky(i) * clrskymax(i)
186 695062 : clrskymax(i) = 1.0_r8
187 : endif
188 493970400 : clrskymax(i) = min(clrskymax(i),1.0_r8-cld(i,k))
189 : endif
190 : end do
191 : end do
192 15507360 : if (ityp == 1) cldtot(1:ncol) = 1.0_r8 - (clrsky(1:ncol) * clrskymax(1:ncol))
193 15507360 : if (ityp == 2) cldlow(1:ncol) = 1.0_r8 - (clrsky(1:ncol) * clrskymax(1:ncol))
194 15507360 : if (ityp == 3) cldmed(1:ncol) = 1.0_r8 - (clrsky(1:ncol) * clrskymax(1:ncol))
195 16256592 : if (ityp == 4) cldhgh(1:ncol) = 1.0_r8 - (clrsky(1:ncol) * clrskymax(1:ncol))
196 : end do
197 :
198 749232 : return
199 : end subroutine cldsav
200 :
201 : end module cloud_cover_diags
|