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

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

Mars GCM:
Cosmetic/practical changes:

  • swmain and lwmain become modules swmain_mod, lwmain_mod
  • Addition of the intent in/out characteristics of variables in swmain_mod and lwmain_mod subroutines

Correction:

  • in callsedim_mod, declaration of variable tau(ngrid,nlay) corrected to tau(ngrid,naerkind)

MV

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