File: C:\NOAA\NEMS_11731\src\atmos\phys\dcyc2.f

1     ! ===================================================================== !
2     !  description:                                                         !
3     !                                                                       !
4     !    dcyc2t3 fits radiative fluxes and heating rates from a coarse      !
5     !    radiation calc time interval into model's more frequent time steps.!
6     !    solar heating rates and fluxes are scaled by the ratio of cosine   !
7     !    of zenith angle at the current time to the mean value used in      !
8     !    radiation calc.  surface downward lw flux is scaled by the ratio   !
9     !    of current surface air temperature (temp**4) to the corresponding  !
10     !    temperature saved during lw radiation calculation. upward lw flux  !
11     !    at the surface is computed by current ground surface temperature.  !
12     !    surface emissivity effect will be taken in other part of the model.!
13     !                                                                       !
14     !  usage:                                                               !
15     !                                                                       !
16     !    call dcyc2t3                                                       !
17     !      inputs:                                                          !
18     !          ( solhr,slag,sdec,cdec,sinlat,coslat,                        !
19     !            xlon,coszen,tsea,tf,tsflw,                                 !
20     !            sfcdsw,sfcnsw,sfcdlw,swh,hlw,                              !
21     !            ix, im, levs,                                              !
22     !      input/output:                                                    !
23     !            dtdt,                                                      !
24     !      outputs:                                                         !
25     !            adjsfcdsw,adjsfcnsw,adjsfcdlw,adjsfculw,xmu,xcosz)         !
26     !                                                                       !
27     !                                                                       !
28     !  program history:                                                     !
29     !          198?  nmc mrf    - created, similar as treatment in gfdl     !
30     !                             radiation treatment                       !
31     !          1994  y. hou     - modified solar zenith angle calculation   !
32     !     nov  2004  x. wu      - add sfc sw downward flux to the variable  !
33     !                             list for sea-ice model                    !
34     !     mar  2008  y. hou     - add cosine of zenith angle as output for  !
35     !                             sunshine duration time calc.              !
36     !     sep  2008  y. hou     - separate net sw and downward lw in slrad, !
37     !                 changed the sign of sfc net sw to consistent with     !
38     !                 other parts of the mdl (positive value defines from   !
39     !                 atmos to the ground). rename output fluxes as adjusted!
40     !                 fluxes. other minor changes such as renaming some of  !
41     !                 passing argument names to be consistent with calling  !
42     !                 program.                                              !
43     !     apr  2009  y. hou     - integrated with the new parallel model    !
44     !                 along with other modifications                        !
45     !                                                                       !
46     !  subprograms called:  none                                            !
47     !                                                                       !
48     !                                                                       !
49     !  ====================  defination of variables  ====================  !
50     !                                                                       !
51     !  inputs:                                                              !
52     !     solhr        - real, forecast time in 24-hour form (hr)           !
53     !     slag         - real, equation of time in radians                  !
54     !     sdec, cdec   - real, sin and cos of the solar declination angle   !
55     !     sinlat(im), coslat(im):                                           !
56     !                  - real, sin and cos of latitude                      !
57     !     xlon   (im)  - real, longitude in radians                         !
58     !     coszen (im)  - real, avg of cosz over daytime sw call interval    !
59     !     tsea   (im)  - real, ground surface temperature (k)               !
60     !     tf     (im)  - real, surface air (layer 1) temperature (k)        !
61     !     tsflw  (im)  - real, sfc air (layer 1) temp in k saved in lw call !
62     !     sfcdsw (im)  - real, total sky sfc downward sw flux ( w/m**2 )    !
63     !     sfcnsw (im)  - real, total sky sfc net sw into ground (w/m**2)    !
64     !     sfcdlw (im)  - real, total sky sfc downward lw flux ( w/m**2 )    !
65     !     swh(ix,levs) - real, total sky sw heating rates ( k/s )           !
66     !     hlw(ix,levs) - real, total sky lw heating rates ( k/s )           !
67     !     ix, im       - integer, horiz. dimention and num of used points   !
68     !     levs         - integer, vertical layer dimension                  !
69     !                                                                       !
70     !  input/output:                                                        !
71     !     dtdt(im,levs)- real, model time step adjusted total radiation     !
72     !                          heating rates ( k/s )                        !
73     !                                                                       !
74     !  outputs:                                                             !
75     !     adjsfcdsw(im)- real, time step adjusted sfc dn sw flux (w/m**2)   !
76     !     adjsfcnsw(im)- real, time step adj sfc net sw into ground (w/m**2)!
77     !     adjsfcdlw(im)- real, time step adjusted sfc dn lw flux (w/m**2)   !
78     !     adjsfculw(im)- real, sfc upward lw flux at current time (w/m**2)  !
79     !     xmu   (im)   - real, time step zenith angle adjust factor for sw  !
80     !     xcosz (im)   - real, cosine of zenith angle at current time step  !
81     !                                                                       !
82     !  ====================    end of description    =====================  !
83     
84     !-----------------------------------
85           subroutine dcyc2t3                                                &
86     !...................................
87     !  ---  inputs:
88          &     ( solhr,slag,sdec,cdec,sinlat,coslat,                        &
89          &       xlon,coszen,tsea,tf,tsflw,                                 &
90          &       sfcdsw,sfcnsw,sfcdlw,swh,hlw,                              &
91          &       ix, im, levs,                                              &
92     !  ---  input/output:
93          &       dtdt,                                                      &
94     !  ---  outputs:
95          &       adjsfcdsw,adjsfcnsw,adjsfcdlw,adjsfculw,xmu,xcosz          &
96          &     )
97     !
98           use machine,         only : kind_phys
99           use physcons,        only : con_pi, con_sbc
100     
101           implicit none
102     !
103     !  ---  constant parameters:
104           real(kind=kind_phys), parameter :: f_eps  = 0.0001
105           real(kind=kind_phys), parameter :: hour12 = 12.0
106     
107     !  ---  inputs:
108           integer, intent(in) :: ix, im, levs
109     
110           real(kind=kind_phys), intent(in) :: solhr, slag, cdec, sdec
111     
112           real(kind=kind_phys), dimension(im), intent(in) ::                &
113          &      sinlat, coslat, xlon, coszen, tsea, tf, tsflw, sfcdlw,      &
114          &      sfcdsw, sfcnsw
115     
116           real(kind=kind_phys), dimension(ix,levs), intent(in) :: swh, hlw
117     
118     !  ---  input/output:
119           real(kind=kind_phys), dimension(im,levs), intent(inout) :: dtdt
120     
121     !  ---  outputs:
122           real(kind=kind_phys), dimension(im), intent(out) ::               &
123          &      adjsfcdsw, adjsfcnsw, adjsfcdlw, adjsfculw, xmu, xcosz
124     
125     !  ---  locals:
126           integer :: i, k
127           real(kind=kind_phys) :: cns, ss, cc, ch, tem1, tem2, xlw(im)
128     !
129     !===> ...  begin here
130     !
131     !  --- ...  compute cosine of solar zenith angle for both hemispheres.
132     
133           cns = con_pi * (solhr - hour12) / hour12 + slag
134     
135           do i = 1, im
136             ss     = sinlat(i) * sdec
137             cc     = coslat(i) * cdec
138             ch     = cc * cos( xlon(i)+cns )
139             xmu(i) = ch + ss
140     
141     !  --- ...  for lw time-step adjustment: compute 4th power of the ratio of
142     !           layer 1 temperature over the value during radiation calculation
143     
144             tem1   = tf(i) / tsflw(i)
145             xlw(i) = tem1 * tem1
146           enddo
147     
148     !  --- ...  normalize by average value over radiation period for daytime.
149     
150           do i = 1, im
151             xcosz(i) = xmu(i)
152             if ( xmu(i) > f_eps .and. coszen(i) > f_eps ) then
153               xmu(i) = xmu(i) / coszen(i)
154             else
155               xmu(i) = 0.0
156             endif
157           enddo
158     
159           do i = 1, im
160     !  --- ...  adjust sfc downward lw flux to account for t changes in layer 1.
161     !           1st line is for original version, 2nd one is for updated version
162             adjsfcdlw(i) = sfcdlw(i) * xlw(i) * xlw(i)
163     
164     !  --- ...  compute sfc upward lw flux from current temp,
165     !      note: sfc emiss effect is not appied at this time
166             tem1         = tsea(i) * tsea(i)
167             adjsfculw(i) =  con_sbc * tem1 * tem1
168     
169     !  --- ...  adjust sfc net and downward sw fluxes for zenith angle changes
170             adjsfcnsw(i) = sfcnsw(i) * xmu(i)
171             adjsfcdsw(i) = sfcdsw(i) * xmu(i)
172           enddo
173     
174     !  --- ...  adjust sw heating rates with zenith angle change and
175     !           add with lw heating to temperature tendency
176     
177           do k = 1, levs
178             do i = 1, im
179               dtdt(i,k) = dtdt(i,k) + swh(i,k)*xmu(i) + hlw(i,k)
180             enddo
181           enddo
182     !
183           return
184     !...................................
185           end subroutine dcyc2t3
186     !-----------------------------------
187     
188