source: trunk/LMDZ.MARS/libf/phymars/callradite_mod.F @ 1974

Last change on this file since 1974 was 1974, checked in by mvals, 6 years ago

Mars GCM:
Integration of the detached dust layer parametrizations (rocket dust storm, slope wind lifting, CW, and dust injection scheme, DB).
Still experimental, default behaviour (rdstorm=.false., dustinjection=0) identical to previous revision.
NB: Updated newstart requires an updated "surface.nc" containing the "hmons" field.
EM+MV

File size: 20.8 KB
Line 
1      MODULE callradite_mod
2
3      IMPLICIT NONE
4
5      CONTAINS
6
7      SUBROUTINE callradite(icount,ngrid,nlayer,nq,zday,ls,pq,albedo,
8     $     emis,mu0,pplev,pplay,pt,tsurf,fract,dist_sol,igout,
9     $     dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw,
10     $     fluxtop_sw,tauref,tau,aerosol,dsodust,tauscaling,
11     $     taucloudtes,rdust,rice,nuice,co2ice,rstormdust,
12     $     totstormfract,clearatm,dsords,
13     $     clearsky,totcloudfrac)
14
15      use aeropacity_mod, only: aeropacity
16      use updatereffrad_mod, only: updatereffrad
17      use dimradmars_mod, only: ndomainsz, nflev, nsun, nir
18      use dimradmars_mod, only: naerkind, name_iaer,
19     &            iaer_dust_conrath,iaer_dust_doubleq,
20     &            iaer_dust_submicron,iaer_h2o_ice,
21     &            iaer_stormdust_doubleq
22      use yomlw_h, only: gcp, nlaylte
23      use comcstfi_h, only: g,cpp
24      use time_phylmdz_mod, only: daysec
25      IMPLICIT NONE
26c=======================================================================
27c   subject:
28c   --------
29c   Subroutine designed to call the main canonic
30c   radiative transfer subroutine "lwmain" et "swmain"
31c   to compute radiative heating and cooling rate and
32c   radiative fluxes to the surface.
33c
34c   These calculations are only valid on the part of the atmosphere
35c   where Local Thermal Equilibrium (NLTE) is verified. In practice
36c   The calculations are only performed for the first "nlaylte"
37c   parameters (nlaylte is calculated by subroutine "nlthermeq"
38c   and stored in module "yomlw_h").
39c
40c   The purpose of this subroutine is to:
41c      1) Make some initial calculation at first call
42c      2) Split the calculation in several sub-grid
43c        ("sub-domain") to save memory and
44c        be able run on a workstation at high resolution
45c        The sub-grid size is defined in dimradmars_mod
46c      3) Compute the 3D scattering parameters depending on the
47c        size distribution of the different tracers (added by JBM)
48c      4) call "lwmain" and "swmain"
49c
50c
51c   authors:   
52c   ------
53c   Francois Forget / Christophe Hourdin / J.-B. Madeleine (2009)
54c
55c
56c   3D scattering scheme user's guide (J.-B. Madeleine)
57c   ---------------------------------
58c
59c   This routine has been modified to take into account 3D, time
60c   dependent scattering properties of the aerosols.
61c---- The look-up tables that contain the scattering parameters
62c   of a given tracer, for different sizes, are read by SUAER.F90.
63c   The names of the corresponding ASCII files have to be set in
64c   this subroutine (file_id variable), and files must be in the
65c   directory specified in datafile_mod. Please make sure that the
66c   ASCII files are correctly written, and that the range
67c   of particle sizes is consistent with what you would expect.
68c---- SUAER.F90 is in charge of reading the ASCII files and averaging
69c   the scattering parameters in each GCM channel, using the three last
70c   equations of Forget et al. 1998 (GRL 25, No.7, p.1105-1108).
71c---- These look-up tables, loaded during the firstcall, are then
72c   constantly used by the subroutine "aeroptproperties.F" to compute,
73c   online, the 3D scattering parameters, based on the size distribution
74c   (reffrad and nueffrad) of the different tracers, in each grid box.
75c   These 3D size distributions are loaded by the "updatereffrad.F"
76c   subroutine. A log-normal distribution is then assumed in
77c   "aeroptproperties.F", along with a Gauss-Legendre integration.
78c---- The optical depth at the visible reference wavelength (set in
79c   SUAER.F90, after the file_id variable) is then computed by
80c   the subroutine "aeropacity.F", by using the size and spatial
81c   distribution of the corresponding tracer. This connection has to
82c   be implemented in "aeropacity.F" when adding a new tracer. To do so,
83c   one can use equation 2 of Forget et al. 1998 (Icarus 131, p.302-316).
84c---- The resulting variables "aerosol", "QVISsQREF3d", "omegaVIS3d" and
85c   "gVIS3d" (same in the infrared) are finally used by lwmain.F and
86c   swmain.F to solve the radiative transfer equation.
87c
88c   changes:
89c   -------
90c
91c   > SRL 7/2000
92c   
93c   This version has been modified to only calculate radiative tendencies
94c   over layers 1..NFLEV (set in dimradmars_mod).  Returns zero for higher
95c   layers, if any.
96c   In other routines, nlayer -> nflev.
97c   Routines affected: lwflux, lwi, lwmain, lwxb, lwxd, lwxn.
98c
99c   > J.-B. Madeleine 10W12
100c   This version uses the variable's splitting, which can be usefull
101c     when performing very high resolution simulation like LES.
102c
103c   ----------
104c   Here, solar band#1 is spectral interval between "long1vis" and "long2vis"
105c   set in dimradmars_mod
106c   Here, solar band#2 is spectral interval between "long2vis" and "long3vis"
107c   set in dimradmars_mod
108c
109c   input:
110c   -----
111c   icount                counter of call to subroutine physic by gcm
112c   ngrid                 number of gridpoint of horizontal grid
113c   nlayer                Number of layer
114c   nq                    Number of tracer
115c   ls                    Solar longitude (Ls) , radian
116c   zday                  Date (time since Ls=0, in martian days)
117c   pq(ngrid,nlayer,nq)   Advected fields
118c
119c   albedo (ngrid,2)      hemispheric surface albedo
120c                         albedo (i,1) : mean albedo for solar band#1
121c                                        (see below)
122c                         albedo (i,2) : mean albedo for solar band#2
123c                                        (see below)
124c   emis                  Thermal IR surface emissivity (no unit)
125c   mu0(ngrid)           cos of solar zenith angle
126c                           (=1 when sun at zenith)
127c   pplay(ngrid,nlayer)    pressure (Pa) in the middle of each layer
128c   pplev(ngrid,nlayer+1)  pressure (Pa) at boundaries of each layer
129c   pt(ngrid,nlayer)       atmospheric temperature in each layer (K)
130c   tsurf(ngrid)           surface temperature (K)
131c   fract(ngrid)         day fraction of the time interval
132c                          =1 during the full day ; =0 during the night
133c   declin                 latitude of subsolar point
134c   dist_sol               sun-Mars distance (AU)
135c   igout                  coordinate of analysed point for debugging
136c   reffrad(ngrid,nlayer,naerkind)  Aerosol effective radius
137c   nueffrad(ngrid,nlayer,naerkind) Aerosol effective variance
138
139c
140c  output:
141c  -------
142c dtlw (ngrid,nlayer)       longwave (IR) heating rate (K/s)
143c dtsw(ngrid,nlayer)        shortwave (Solar) heating rate (K/s)
144c fluxsurf_lw(ngrid)        surface downward flux tota LW (thermal IR) (W.m-2)
145c fluxsurf_sw(ngrid,1)      surface downward flux SW for solar band#1 (W.m-2)
146c fluxsurf_sw(ngrid,2)      surface downward flux SW for solar band#2 (W.m-2)
147c
148c fluxtop_lw(ngrid)         outgoing upward flux tota LW (thermal IR) (W.m-2)
149c fluxtop_sw(ngrid,1)       outgoing upward flux SW for solar band#1 (W.m-2)
150c fluxtop_sw(ngrid,2)       outgoing upward flux SW for solar band#2 (W.m-2)
151
152c   tauref       Prescribed mean column optical depth at 610 Pa
153c   tau          Column total visible dust optical depth at each point
154c   aerosol(ngrid,nlayer,naerkind)    aerosol extinction optical depth
155c                         at reference wavelength "longrefvis" set
156c                         in dimradmars_h , in each layer, for one of
157c                         the "naerkind" kind of aerosol optical
158c                         properties.
159c=======================================================================
160c
161c    Declarations :
162c    -------------
163c
164      include "callkeys.h"
165
166c-----------------------------------------------------------------------
167c    Input/Output
168c    ------------
169      INTEGER,INTENT(IN) :: icount       
170      INTEGER,INTENT(IN) :: ngrid,nlayer,nq
171      INTEGER,INTENT(IN) :: igout
172
173      REAL,INTENT(IN) :: pq(ngrid,nlayer,nq)
174      REAL,INTENT(INOUT) :: tauscaling(ngrid) ! Conversion factor for
175                               ! qdust and Ndust
176      REAL,INTENT(IN) :: albedo(ngrid,2),emis(ngrid)
177      REAL,INTENT(IN) :: ls,zday
178
179      REAL,INTENT(IN) :: pplev(ngrid,nlayer+1),pplay(ngrid,nlayer)
180      REAL,INTENT(IN) :: pt(ngrid,nlayer)
181      REAL,INTENT(IN) :: tsurf(ngrid)
182      REAL,INTENT(IN) :: dist_sol,mu0(ngrid),fract(ngrid)
183      REAL,INTENT(OUT) :: dtlw(ngrid,nlayer),dtsw(ngrid,nlayer)
184      REAL,INTENT(OUT) :: fluxsurf_lw(ngrid), fluxtop_lw(ngrid)
185      REAL,INTENT(OUT) :: fluxsurf_sw(ngrid,2), fluxtop_sw(ngrid,2)
186
187      REAL,INTENT(OUT) :: tauref(ngrid), tau(ngrid,naerkind)
188      REAL,INTENT(OUT) :: taucloudtes(ngrid)! Cloud opacity at infrared
189                               !   reference wavelength using
190                               !   Qabs instead of Qext
191                               !   (direct comparison with TES)
192      REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind)
193      REAL,INTENT(OUT) :: dsodust(ngrid,nlayer)
194      REAL,INTENT(OUT) :: rdust(ngrid,nlayer)  ! Dust geometric mean radius (m)
195      REAL,INTENT(OUT) :: rice(ngrid,nlayer)   ! Ice geometric mean radius (m)
196      REAL,INTENT(OUT) :: nuice(ngrid,nlayer)  ! Estimated effective variance
197      REAL,INTENT(IN) :: co2ice(ngrid)           ! co2 ice surface layer (kg.m-2)
198
199c     rocket dust storm
200      LOGICAL,INTENT(IN) :: clearatm ! true for background dust
201      REAL,INTENT(IN) :: totstormfract(ngrid) ! dust storm mesh fraction
202      REAL,INTENT(OUT) :: rstormdust(ngrid,nlayer)  ! Storm dust geometric mean radius (m)
203      REAL dsords(ngrid,nlayer) ! density scaled opacity for rocket dust storm dust
204
205c     sub-grid scale water ice clouds
206      LOGICAL,INTENT(IN) :: clearsky
207      REAL,INTENT(IN) :: totcloudfrac(ngrid)
208
209c
210c    Local variables :
211c    -----------------
212
213      INTEGER j,l,ig,n,ich
214      INTEGER aer_count,iaer
215      INTEGER jd,ig0,nd
216
217      real  cste_mars ! solar constant on Mars (Wm-2)
218      REAL ptlev(ngrid,nlayer+1)
219
220      INTEGER :: ndomain
221
222c     Thermal IR net radiative budget (W m-2)
223
224      real znetrad(ndomainsz,nflev)
225
226      real zfluxd_sw(ndomainsz,nflev+1,2)
227      real zfluxu_sw(ndomainsz,nflev+1,2)
228
229      REAL zplev(ndomainsz,nflev+1)
230      REAL zztlev(ndomainsz,nflev+1)
231      REAL zplay(ndomainsz,nflev)
232      REAL zt(ndomainsz,nflev)
233      REAL zaerosol(ndomainsz,nflev,naerkind)
234      REAL zalbedo(ndomainsz,2)
235      REAL zdp(ndomainsz,nflev)
236      REAL zdt0(ndomainsz)
237
238      REAL zzdtlw(ndomainsz,nflev)
239      REAL zzdtsw(ndomainsz,nflev)
240      REAL zzflux(ndomainsz,6)
241      real zrmuz
242
243      REAL :: zQVISsQREF3d(ndomainsz,nflev,nsun,naerkind)
244      REAL :: zomegaVIS3d(ndomainsz,nflev,nsun,naerkind)
245      REAL :: zgVIS3d(ndomainsz,nflev,nsun,naerkind)
246
247      REAL :: zQIRsQREF3d(ndomainsz,nflev,nir,naerkind)
248      REAL :: zomegaIR3d(ndomainsz,nflev,nir,naerkind)
249      REAL :: zgIR3d(ndomainsz,nflev,nir,naerkind)
250
251c     Aerosol size distribution
252      REAL :: reffrad(ngrid,nlayer,naerkind)
253      REAL :: nueffrad(ngrid,nlayer,naerkind)
254c     Aerosol optical properties
255      REAL :: QVISsQREF3d(ngrid,nlayer,nsun,naerkind)
256      REAL :: omegaVIS3d(ngrid,nlayer,nsun,naerkind)
257      REAL :: gVIS3d(ngrid,nlayer,nsun,naerkind)
258
259      REAL :: QIRsQREF3d(ngrid,nlayer,nir,naerkind)
260      REAL :: omegaIR3d(ngrid,nlayer,nir,naerkind)
261      REAL :: gIR3d(ngrid,nlayer,nir,naerkind)
262
263      REAL :: QREFvis3d(ngrid,nlayer,naerkind)
264      REAL :: QREFir3d(ngrid,nlayer,naerkind)
265
266      REAL :: omegaREFvis3d(ngrid,nlayer,naerkind)
267      REAL :: omegaREFir3d(ngrid,nlayer,naerkind)
268
269c   local saved variables
270c   ---------------------
271
272      real zco2   ! volume fraction of CO2 in Mars atmosphere
273      DATA zco2/0.95/
274      SAVE zco2
275
276      LOGICAL firstcall
277      DATA firstcall/.true./
278      SAVE firstcall
279
280c----------------------------------------------------------------------
281
282c     Initialisation
283c     --------------
284
285! compute ndomain
286! AS: moved out of firstcall to allow nesting+evoluting domain
287! ------------------------------------------------------------
288      ndomain= (ngrid-1) / ndomainsz + 1
289
290      IF (firstcall) THEN
291
292         write(*,*) 'Splitting radiative calculations: ',
293     $              ' ngrid,ndomainsz,ndomain',
294     $                ngrid,ndomainsz,ndomain
295
296c        Assign a number to the different scatterers
297c        -------------------------------------------
298
299         iaer_dust_conrath=0
300         iaer_dust_doubleq=0
301         iaer_dust_submicron=0
302         iaer_h2o_ice=0
303         iaer_stormdust_doubleq=0
304
305         aer_count=0
306         if (.NOT.active) then
307           do iaer=1,naerkind
308             if (name_iaer(iaer).eq."dust_conrath") then
309               iaer_dust_conrath = iaer
310               aer_count = aer_count + 1
311             endif
312           enddo
313         endif
314         if (doubleq.AND.active) then
315           do iaer=1,naerkind
316             if (name_iaer(iaer).eq."dust_doubleq") then
317               iaer_dust_doubleq = iaer
318               aer_count = aer_count + 1
319             endif
320           enddo
321         endif
322         if (submicron.AND.active) then
323           do iaer=1,naerkind
324             if (name_iaer(iaer).eq."dust_submicron") then
325               iaer_dust_submicron = iaer
326               aer_count = aer_count + 1
327             endif
328           enddo
329         endif
330         if (water.AND.activice) then
331           do iaer=1,naerkind
332             if (name_iaer(iaer).eq."h2o_ice") then
333               iaer_h2o_ice = iaer
334               aer_count = aer_count + 1
335             endif
336           enddo
337         endif
338         if (rdstorm.AND.active) then
339           do iaer=1,naerkind
340             if (name_iaer(iaer).eq."stormdust_doubleq") then
341               iaer_stormdust_doubleq = iaer
342               aer_count = aer_count + 1
343             endif
344           enddo
345         end if
346
347c        Check that we identified all tracers:
348         if (aer_count.ne.naerkind) then
349           write(*,*) "callradite: found only ",aer_count," scatterers"
350           write(*,*) "               expected ",naerkind
351           write(*,*) "please make sure that the number of"
352           write(*,*) "scatterers in scatterers.h, the names"
353           write(*,*) "in callradite.F, and the flags in"
354           write(*,*) "callphys.def are all consistent!"
355           do iaer=1,naerkind
356             write(*,*)'      ',iaer,' ',trim(name_iaer(iaer))
357           enddo
358           stop
359         else
360           write(*,*) "callradite: found all scatterers, namely:"
361           do iaer=1,naerkind
362             write(*,*)'      ',iaer,' ',trim(name_iaer(iaer))
363           enddo
364         endif
365c        -------------------------------------------
366
367         gcp = g/cpp
368
369c        Loading the optical properties in external look-up tables:
370         CALL SUAER
371!         CALL SULW ! this step is now done in ini_yomlw_h
372
373         if (ngrid .EQ. 1) then
374           if (ndomainsz .NE. 1) then
375             print*
376             print*,'ATTENTION !!!'
377             print*,'pour tourner en 1D, '
378             print*,'fixer ndomainsz=1 dans phymars/dimradmars_h'
379             print*
380             call exit(1)
381           endif
382         endif
383
384         firstcall=.false.
385      END IF
386
387c     Computing aerosol optical properties and opacity
388c     ------------------------------------------------
389
390c     Updating aerosol size distributions:
391      CALL updatereffrad(ngrid,nlayer,
392     &                rdust,rstormdust,rice,nuice,
393     &                reffrad,nueffrad,
394     &                pq,tauscaling,tau,pplay)
395
396c     Computing 3D scattering parameters:
397      CALL aeroptproperties(ngrid,nlayer,reffrad,nueffrad,
398     &                      QVISsQREF3d,omegaVIS3d,gVIS3d,
399     &                      QIRsQREF3d,omegaIR3d,gIR3d,
400     &                      QREFvis3d,QREFir3d,
401     &                      omegaREFvis3d,omegaREFir3d)
402
403c     Computing aerosol optical depth in each layer:
404      CALL aeropacity(ngrid,nlayer,nq,zday,pplay,pplev,ls,
405     &    pq,tauscaling,tauref,tau,taucloudtes,aerosol,dsodust,reffrad,
406     &    QREFvis3d,QREFir3d,omegaREFir3d,
407     &    totstormfract,clearatm,dsords,
408     &    clearsky,totcloudfrac)
409
410c     Starting loop on sub-domain
411c     ----------------------------
412
413      DO jd=1,ndomain
414        ig0=(jd-1)*ndomainsz
415        if (jd.eq.ndomain) then
416         nd=ngrid-ig0
417        else
418         nd=ndomainsz
419        endif
420
421c       Spliting input variable in sub-domain input variables
422c       ---------------------------------------------------
423
424        do l=1,nlaylte
425         do ig = 1,nd
426           do iaer = 1, naerkind
427             do ich = 1, nsun
428               zQVISsQREF3d(ig,l,ich,iaer) =
429     &                           QVISsQREF3d(ig0+ig,l,ich,iaer)
430               zomegaVIS3d(ig,l,ich,iaer) =
431     &                           omegaVIS3d(ig0+ig,l,ich,iaer)
432               zgVIS3d(ig,l,ich,iaer) =
433     &                           gVIS3d(ig0+ig,l,ich,iaer)
434             enddo
435             do ich = 1, nir
436               zQIRsQREF3d(ig,l,ich,iaer) =
437     &                           QIRsQREF3d(ig0+ig,l,ich,iaer)
438               zomegaIR3d(ig,l,ich,iaer) =
439     &                           omegaIR3d(ig0+ig,l,ich,iaer)
440               zgIR3d(ig,l,ich,iaer) =
441     &                           gIR3d(ig0+ig,l,ich,iaer)
442             enddo
443           enddo
444         enddo
445        enddo
446
447        do l=1,nlaylte+1
448         do ig = 1,nd
449          zplev(ig,l) = pplev(ig0+ig,l)
450         enddo
451        enddo
452
453        do l=1,nlaylte
454         do ig = 1,nd
455          zplay(ig,l) = pplay(ig0+ig,l)
456          zt(ig,l) = pt(ig0+ig,l)
457c         Thickness of each layer (Pa) :
458          zdp(ig,l)= pplev(ig0+ig,l) - pplev(ig0+ig,l+1)
459         enddo
460        enddo
461
462        do n=1,naerkind
463          do l=1,nlaylte
464            do ig=1,nd
465              zaerosol(ig,l,n) = aerosol(ig0+ig,l,n)
466            enddo
467          enddo
468        enddo
469
470        do j=1,2
471          do ig = 1,nd
472           zalbedo(ig,j) = albedo(ig0+ig,j)
473          enddo
474        enddo
475
476c       Intermediate  levels: (computing tlev)
477c       ---------------------------------------
478c       Extrapolation for the air temperature above the surface
479        DO ig=1,nd
480              zztlev(ig,1)=zt(ig,1)+
481     s        (zplev(ig,1)-zplay(ig,1))*
482     s        (zt(ig,1)-zt(ig,2))/(zplay(ig,1)-zplay(ig,2))
483
484              zdt0(ig) = tsurf(ig0+ig) - zztlev(ig,1)
485        ENDDO
486
487        DO l=2,nlaylte
488         DO ig=1, nd
489               zztlev(ig,l)=0.5*(zt(ig,l-1)+zt(ig,l))
490         ENDDO
491        ENDDO
492
493        DO ig=1, nd
494           zztlev(ig,nlaylte+1)=zt(ig,nlaylte)
495        ENDDO
496
497
498c       Longwave ("lw") radiative transfer (= thermal infrared)
499c       -------------------------------------------------------
500        call lwmain (ig0,icount,nd,nflev
501     .        ,zdp,zdt0,emis(ig0+1),zplev,zztlev,zt
502     .        ,zaerosol,zzdtlw
503     .        ,fluxsurf_lw(ig0+1),fluxtop_lw(ig0+1)
504     .        ,znetrad
505     &        ,zQIRsQREF3d,zomegaIR3d,zgIR3d
506     &        ,co2ice(ig0+1))
507
508c       Shortwave ("sw") radiative transfer (= solar radiation)
509c       -------------------------------------------------------
510c          Mars solar constant (W m-2)
511c          1370 W.m-2 is the solar constant at 1 AU.
512           cste_mars=1370./(dist_sol*dist_sol)
513
514           call swmain ( nd, nflev,
515     S     cste_mars, zalbedo,
516     S     mu0(ig0+1), zdp, zplev, zaerosol, fract(ig0+1),
517     S     zzdtsw, zfluxd_sw, zfluxu_sw,
518     &     zQVISsQREF3d,zomegaVIS3d,zgVIS3d)
519
520c       ------------------------------------------------------------
521c       Un-spliting output variable from sub-domain input variables
522c       ------------------------------------------------------------
523
524        do l=1,nlaylte
525         do ig = 1,nd
526          dtlw(ig0+ig,l) = zzdtlw(ig,l)
527          dtsw(ig0+ig,l) = zzdtsw(ig,l)
528         enddo
529        enddo
530
531        do l=1,nlaylte+1
532         do ig = 1,nd
533          ptlev(ig0+ig,l) = zztlev(ig,l)
534         enddo
535        enddo
536
537        do ig = 1,nd
538          fluxsurf_sw(ig0+ig,1) = zfluxd_sw(ig,1,1)
539          fluxsurf_sw(ig0+ig,2) = zfluxd_sw(ig,1,2)
540          fluxtop_sw(ig0+ig,1) = zfluxu_sw(ig,nlaylte+1,1)
541          fluxtop_sw(ig0+ig,2) = zfluxu_sw(ig,nlaylte+1,2)
542        enddo
543
544      ENDDO         !   (boucle jd=1, ndomain)
545
546c     Zero tendencies for any remaining layers between nlaylte and nlayer
547      if (nlayer.gt.nlaylte) then
548         do l = nlaylte+1, nlayer
549            do ig = 1, ngrid
550               dtlw(ig, l) = 0.
551               dtsw(ig, l) = 0.
552            enddo
553         enddo
554      endif
555
556c     Output for debugging if lwrite=T
557c     --------------------------------
558c     Write all nlayer layers, even though only nlaylte layers may have
559c     non-zero tendencies.
560
561         IF(lwrite) THEN
562            PRINT*,'Diagnotique for the radiation'
563            PRINT*,'albedo, emissiv, mu0,fract,fluxsurf_lw,fluxsurf_sw'
564            PRINT*,albedo(igout,1),emis(igout),mu0(igout),
565     s           fract(igout), fluxsurf_lw(igout),
566     $     fluxsurf_sw(igout,1)+fluxsurf_sw(igout,2)
567            PRINT*,'Tlay Tlev Play Plev dT/dt SW dT/dt LW (K/s)'
568            PRINT*,'daysec',daysec
569            DO l=1,nlayer
570               PRINT*,pt(igout,l),ptlev(igout,l),
571     s         pplay(igout,l),pplev(igout,l),
572     s         dtsw(igout,l),dtlw(igout,l)
573            ENDDO
574         ENDIF
575
576
577      END SUBROUTINE callradite
578
579      END MODULE callradite_mod
Note: See TracBrowser for help on using the repository browser.