source: trunk/mars/libf/phymars/callradite.F @ 132

Last change on this file since 132 was 131, checked in by aslmd, 14 years ago

mars

M 130 mars/libf/phymars/callradite.F
M 130 mars/README
no more need to modify callradite.F prior to compilation
[but still dimradmars.h must be modified]

--> in callradite.F we now have

-- DEFAULT name_iaer(1) is "dust_conrath"
-- IF (doubleq.AND.active) name_iaer(1) = "dust_doubleq"
-- IF (water.AND.activice) name_iaer(2) = "h2o_ice"

M 130 000-USERS
small meeting about SVN for "planeto" team

M 130 mars/libf/phymars/meso_physiq.F
purely cosmetic changes

M 130 mesoscale/PLOT/MINIMAL/map_latlon.pro
M 130 mesoscale/PLOT/SPEC/GW/gravitwave2.pro
M 130 mesoscale/PLOT/SPEC/GW/gravitwaveprof.pro
M 130 mesoscale/PLOT/SPEC/GW/gravitwave_inc.pro
A 0 mesoscale/PLOT/SPEC/MAP/defs/THARSIS_newphys_bugnoRAC_RACgcm_TAUTES.map_uvt_inc.pro
M 130 mesoscale/PLOT/SPEC/MAP/map_uvt.pro
graphic stuff for GW studies and H2O cycle

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