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

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

Mars GCM:
Implementation of a new parametrization of the dust entrainment by slope winds above the sub-grid scale topography. The parametrization is activated with the flag slpwind=.true. (set to "false" by
default) in callphys.def. The new parametrization involves the new tracers topdust_mass and topdust_number.
MV

File size: 21.5 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,rtopdust,
12     $     totstormfract,clearatm,dsords,alpha_hmons,nohmons,
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,iaer_topdust_doubleq
22      use yomlw_h, only: gcp, nlaylte
23      use comcstfi_h, only: g,cpp
24      use time_phylmdz_mod, only: daysec
25      use lwmain_mod, only: lwmain
26      use swmain_mod, only: swmain
27      IMPLICIT NONE
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"
40c   and stored in module "yomlw_h").
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
47c        The sub-grid size is defined in dimradmars_mod
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
67c   directory specified in datafile_mod. Please make sure that the
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
96c   over layers 1..NFLEV (set in dimradmars_mod).  Returns zero for higher
97c   layers, if any.
98c   In other routines, nlayer -> nflev.
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"
107c   set in dimradmars_mod
108c   Here, solar band#2 is spectral interval between "long2vis" and "long3vis"
109c   set in dimradmars_mod
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)
127c   mu0(ngrid)           cos of solar zenith angle
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)
133c   fract(ngrid)         day fraction of the time interval
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
154c   tauref       Prescribed mean column optical depth at 610 Pa
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
158c                         in dimradmars_h , in each layer, for one of
159c                         the "naerkind" kind of aerosol optical
160c                         properties.
161c=======================================================================
162c
163c    Declarations :
164c    -------------
165c
166      include "callkeys.h"
167
168c-----------------------------------------------------------------------
169c    Input/Output
170c    ------------
171      INTEGER,INTENT(IN) :: icount       
172      INTEGER,INTENT(IN) :: ngrid,nlayer,nq
173      INTEGER,INTENT(IN) :: igout
174
175      REAL,INTENT(IN) :: pq(ngrid,nlayer,nq)
176      REAL,INTENT(INOUT) :: tauscaling(ngrid) ! Conversion factor for
177                               ! qdust and Ndust
178      REAL,INTENT(IN) :: albedo(ngrid,2),emis(ngrid)
179      REAL,INTENT(IN) :: ls,zday
180
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)
188
189      REAL,INTENT(OUT) :: tauref(ngrid), tau(ngrid,naerkind)
190      REAL,INTENT(OUT) :: taucloudtes(ngrid)! Cloud opacity at infrared
191                               !   reference wavelength using
192                               !   Qabs instead of Qext
193                               !   (direct comparison with TES)
194      REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind)
195      REAL,INTENT(OUT) :: dsodust(ngrid,nlayer)
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)
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
207c     entrainment by slope wind
208      LOGICAL, INTENT(IN) :: nohmons ! true for background dust
209      REAL, INTENT(IN) :: alpha_hmons(ngrid) ! sub-grid scale topography mesh fraction
210      REAL,INTENT(OUT) :: rtopdust(ngrid,nlayer)  ! Topdust geometric mean radius (m)
211
212c     sub-grid scale water ice clouds
213      LOGICAL,INTENT(IN) :: clearsky
214      REAL,INTENT(IN) :: totcloudfrac(ngrid)
215
216c
217c    Local variables :
218c    -----------------
219
220      INTEGER j,l,ig,n,ich
221      INTEGER aer_count,iaer
222      INTEGER jd,ig0,nd
223
224      real  cste_mars ! solar constant on Mars (Wm-2)
225      REAL ptlev(ngrid,nlayer+1)
226
227      INTEGER :: ndomain
228
229c     Thermal IR net radiative budget (W m-2)
230
231      real znetrad(ndomainsz,nflev)
232
233      real zfluxd_sw(ndomainsz,nflev+1,2)
234      real zfluxu_sw(ndomainsz,nflev+1,2)
235
236      REAL zplev(ndomainsz,nflev+1)
237      REAL zztlev(ndomainsz,nflev+1)
238      REAL zplay(ndomainsz,nflev)
239      REAL zt(ndomainsz,nflev)
240      REAL zaerosol(ndomainsz,nflev,naerkind)
241      REAL zalbedo(ndomainsz,2)
242      REAL zdp(ndomainsz,nflev)
243      REAL zdt0(ndomainsz)
244
245      REAL zzdtlw(ndomainsz,nflev)
246      REAL zzdtsw(ndomainsz,nflev)
247      REAL zzflux(ndomainsz,6)
248      real zrmuz
249
250      REAL :: zQVISsQREF3d(ndomainsz,nflev,nsun,naerkind)
251      REAL :: zomegaVIS3d(ndomainsz,nflev,nsun,naerkind)
252      REAL :: zgVIS3d(ndomainsz,nflev,nsun,naerkind)
253
254      REAL :: zQIRsQREF3d(ndomainsz,nflev,nir,naerkind)
255      REAL :: zomegaIR3d(ndomainsz,nflev,nir,naerkind)
256      REAL :: zgIR3d(ndomainsz,nflev,nir,naerkind)
257
258c     Aerosol size distribution
259      REAL :: reffrad(ngrid,nlayer,naerkind)
260      REAL :: nueffrad(ngrid,nlayer,naerkind)
261c     Aerosol optical properties
262      REAL :: QVISsQREF3d(ngrid,nlayer,nsun,naerkind)
263      REAL :: omegaVIS3d(ngrid,nlayer,nsun,naerkind)
264      REAL :: gVIS3d(ngrid,nlayer,nsun,naerkind)
265
266      REAL :: QIRsQREF3d(ngrid,nlayer,nir,naerkind)
267      REAL :: omegaIR3d(ngrid,nlayer,nir,naerkind)
268      REAL :: gIR3d(ngrid,nlayer,nir,naerkind)
269
270      REAL :: QREFvis3d(ngrid,nlayer,naerkind)
271      REAL :: QREFir3d(ngrid,nlayer,naerkind)
272
273      REAL :: omegaREFvis3d(ngrid,nlayer,naerkind)
274      REAL :: omegaREFir3d(ngrid,nlayer,naerkind)
275
276c   local saved variables
277c   ---------------------
278
279      real zco2   ! volume fraction of CO2 in Mars atmosphere
280      DATA zco2/0.95/
281      SAVE zco2
282
283      LOGICAL firstcall
284      DATA firstcall/.true./
285      SAVE firstcall
286
287c----------------------------------------------------------------------
288
289c     Initialisation
290c     --------------
291
292! compute ndomain
293! AS: moved out of firstcall to allow nesting+evoluting domain
294! ------------------------------------------------------------
295      ndomain= (ngrid-1) / ndomainsz + 1
296
297      IF (firstcall) THEN
298
299         write(*,*) 'Splitting radiative calculations: ',
300     $              ' ngrid,ndomainsz,ndomain',
301     $                ngrid,ndomainsz,ndomain
302
303c        Assign a number to the different scatterers
304c        -------------------------------------------
305
306         iaer_dust_conrath=0
307         iaer_dust_doubleq=0
308         iaer_dust_submicron=0
309         iaer_h2o_ice=0
310         iaer_stormdust_doubleq=0
311         iaer_topdust_doubleq=0
312
313         aer_count=0
314         if (.NOT.active) then
315           do iaer=1,naerkind
316             if (name_iaer(iaer).eq."dust_conrath") then
317               iaer_dust_conrath = iaer
318               aer_count = aer_count + 1
319             endif
320           enddo
321         endif
322         if (doubleq.AND.active) then
323           do iaer=1,naerkind
324             if (name_iaer(iaer).eq."dust_doubleq") then
325               iaer_dust_doubleq = iaer
326               aer_count = aer_count + 1
327             endif
328           enddo
329         endif
330         if (submicron.AND.active) then
331           do iaer=1,naerkind
332             if (name_iaer(iaer).eq."dust_submicron") then
333               iaer_dust_submicron = iaer
334               aer_count = aer_count + 1
335             endif
336           enddo
337         endif
338         if (water.AND.activice) then
339           do iaer=1,naerkind
340             if (name_iaer(iaer).eq."h2o_ice") then
341               iaer_h2o_ice = iaer
342               aer_count = aer_count + 1
343             endif
344           enddo
345         endif
346         if (rdstorm.AND.active) then
347           do iaer=1,naerkind
348             if (name_iaer(iaer).eq."stormdust_doubleq") then
349               iaer_stormdust_doubleq = iaer
350               aer_count = aer_count + 1
351             endif
352           enddo
353         end if
354         if (slpwind.AND.active) then
355           do iaer=1,naerkind
356             if (name_iaer(iaer).eq."topdust_doubleq") then
357               iaer_topdust_doubleq = iaer
358               aer_count = aer_count + 1
359             endif
360           enddo
361         end if
362
363c        Check that we identified all tracers:
364         if (aer_count.ne.naerkind) then
365           write(*,*) "callradite: found only ",aer_count," scatterers"
366           write(*,*) "               expected ",naerkind
367           write(*,*) "please make sure that the number of"
368           write(*,*) "scatterers in scatterers.h, the names"
369           write(*,*) "in callradite.F, and the flags in"
370           write(*,*) "callphys.def are all consistent!"
371           do iaer=1,naerkind
372             write(*,*)'      ',iaer,' ',trim(name_iaer(iaer))
373           enddo
374           stop
375         else
376           write(*,*) "callradite: found all scatterers, namely:"
377           do iaer=1,naerkind
378             write(*,*)'      ',iaer,' ',trim(name_iaer(iaer))
379           enddo
380         endif
381c        -------------------------------------------
382
383         gcp = g/cpp
384
385c        Loading the optical properties in external look-up tables:
386         CALL SUAER
387!         CALL SULW ! this step is now done in ini_yomlw_h
388
389         if (ngrid .EQ. 1) then
390           if (ndomainsz .NE. 1) then
391             print*
392             print*,'ATTENTION !!!'
393             print*,'pour tourner en 1D, '
394             print*,'fixer ndomainsz=1 dans phymars/dimradmars_h'
395             print*
396             call exit(1)
397           endif
398         endif
399
400         firstcall=.false.
401      END IF
402
403c     Computing aerosol optical properties and opacity
404c     ------------------------------------------------
405
406c     Updating aerosol size distributions:
407      CALL updatereffrad(ngrid,nlayer,
408     &                rdust,rstormdust,rtopdust,rice,nuice,
409     &                reffrad,nueffrad,
410     &                pq,tauscaling,tau,pplay)
411
412c     Computing 3D scattering parameters:
413      CALL aeroptproperties(ngrid,nlayer,reffrad,nueffrad,
414     &                      QVISsQREF3d,omegaVIS3d,gVIS3d,
415     &                      QIRsQREF3d,omegaIR3d,gIR3d,
416     &                      QREFvis3d,QREFir3d,
417     &                      omegaREFvis3d,omegaREFir3d)
418
419c     Computing aerosol optical depth in each layer:
420      CALL aeropacity(ngrid,nlayer,nq,zday,pplay,pplev,ls,
421     &    pq,tauscaling,tauref,tau,taucloudtes,aerosol,dsodust,reffrad,
422     &    QREFvis3d,QREFir3d,omegaREFir3d,
423     &    totstormfract,clearatm,dsords,
424     &    alpha_hmons,nohmons,
425     &    clearsky,totcloudfrac)
426
427c     Starting loop on sub-domain
428c     ----------------------------
429
430      DO jd=1,ndomain
431        ig0=(jd-1)*ndomainsz
432        if (jd.eq.ndomain) then
433         nd=ngrid-ig0
434        else
435         nd=ndomainsz
436        endif
437
438c       Spliting input variable in sub-domain input variables
439c       ---------------------------------------------------
440
441        do l=1,nlaylte
442         do ig = 1,nd
443           do iaer = 1, naerkind
444             do ich = 1, nsun
445               zQVISsQREF3d(ig,l,ich,iaer) =
446     &                           QVISsQREF3d(ig0+ig,l,ich,iaer)
447               zomegaVIS3d(ig,l,ich,iaer) =
448     &                           omegaVIS3d(ig0+ig,l,ich,iaer)
449               zgVIS3d(ig,l,ich,iaer) =
450     &                           gVIS3d(ig0+ig,l,ich,iaer)
451             enddo
452             do ich = 1, nir
453               zQIRsQREF3d(ig,l,ich,iaer) =
454     &                           QIRsQREF3d(ig0+ig,l,ich,iaer)
455               zomegaIR3d(ig,l,ich,iaer) =
456     &                           omegaIR3d(ig0+ig,l,ich,iaer)
457               zgIR3d(ig,l,ich,iaer) =
458     &                           gIR3d(ig0+ig,l,ich,iaer)
459             enddo
460           enddo
461         enddo
462        enddo
463
464        do l=1,nlaylte+1
465         do ig = 1,nd
466          zplev(ig,l) = pplev(ig0+ig,l)
467         enddo
468        enddo
469
470        do l=1,nlaylte
471         do ig = 1,nd
472          zplay(ig,l) = pplay(ig0+ig,l)
473          zt(ig,l) = pt(ig0+ig,l)
474c         Thickness of each layer (Pa) :
475          zdp(ig,l)= pplev(ig0+ig,l) - pplev(ig0+ig,l+1)
476         enddo
477        enddo
478
479        do n=1,naerkind
480          do l=1,nlaylte
481            do ig=1,nd
482              zaerosol(ig,l,n) = aerosol(ig0+ig,l,n)
483            enddo
484          enddo
485        enddo
486
487        do j=1,2
488          do ig = 1,nd
489           zalbedo(ig,j) = albedo(ig0+ig,j)
490          enddo
491        enddo
492
493c       Intermediate  levels: (computing tlev)
494c       ---------------------------------------
495c       Extrapolation for the air temperature above the surface
496        DO ig=1,nd
497              zztlev(ig,1)=zt(ig,1)+
498     s        (zplev(ig,1)-zplay(ig,1))*
499     s        (zt(ig,1)-zt(ig,2))/(zplay(ig,1)-zplay(ig,2))
500
501              zdt0(ig) = tsurf(ig0+ig) - zztlev(ig,1)
502        ENDDO
503
504        DO l=2,nlaylte
505         DO ig=1, nd
506               zztlev(ig,l)=0.5*(zt(ig,l-1)+zt(ig,l))
507         ENDDO
508        ENDDO
509
510        DO ig=1, nd
511           zztlev(ig,nlaylte+1)=zt(ig,nlaylte)
512        ENDDO
513
514
515c       Longwave ("lw") radiative transfer (= thermal infrared)
516c       -------------------------------------------------------
517        call lwmain (ig0,icount,nd,nflev
518     .        ,zdp,zdt0,emis(ig0+1),zplev,zztlev,zt
519     .        ,zaerosol,zzdtlw
520     .        ,fluxsurf_lw(ig0+1),fluxtop_lw(ig0+1)
521     .        ,znetrad
522     &        ,zQIRsQREF3d,zomegaIR3d,zgIR3d
523     &        ,co2ice(ig0+1))
524
525c       Shortwave ("sw") radiative transfer (= solar radiation)
526c       -------------------------------------------------------
527c          Mars solar constant (W m-2)
528c          1370 W.m-2 is the solar constant at 1 AU.
529           cste_mars=1370./(dist_sol*dist_sol)
530
531           call swmain ( nd, nflev,
532     S     cste_mars, zalbedo,
533     S     mu0(ig0+1), zdp, zplev, zaerosol, fract(ig0+1),
534     S     zzdtsw, zfluxd_sw, zfluxu_sw,
535     &     zQVISsQREF3d,zomegaVIS3d,zgVIS3d)
536
537c       ------------------------------------------------------------
538c       Un-spliting output variable from sub-domain input variables
539c       ------------------------------------------------------------
540
541        do l=1,nlaylte
542         do ig = 1,nd
543          dtlw(ig0+ig,l) = zzdtlw(ig,l)
544          dtsw(ig0+ig,l) = zzdtsw(ig,l)
545         enddo
546        enddo
547
548        do l=1,nlaylte+1
549         do ig = 1,nd
550          ptlev(ig0+ig,l) = zztlev(ig,l)
551         enddo
552        enddo
553
554        do ig = 1,nd
555          fluxsurf_sw(ig0+ig,1) = zfluxd_sw(ig,1,1)
556          fluxsurf_sw(ig0+ig,2) = zfluxd_sw(ig,1,2)
557          fluxtop_sw(ig0+ig,1) = zfluxu_sw(ig,nlaylte+1,1)
558          fluxtop_sw(ig0+ig,2) = zfluxu_sw(ig,nlaylte+1,2)
559        enddo
560
561      ENDDO         !   (boucle jd=1, ndomain)
562
563c     Zero tendencies for any remaining layers between nlaylte and nlayer
564      if (nlayer.gt.nlaylte) then
565         do l = nlaylte+1, nlayer
566            do ig = 1, ngrid
567               dtlw(ig, l) = 0.
568               dtsw(ig, l) = 0.
569            enddo
570         enddo
571      endif
572
573c     Output for debugging if lwrite=T
574c     --------------------------------
575c     Write all nlayer layers, even though only nlaylte layers may have
576c     non-zero tendencies.
577
578         IF(lwrite) THEN
579            PRINT*,'Diagnotique for the radiation'
580            PRINT*,'albedo, emissiv, mu0,fract,fluxsurf_lw,fluxsurf_sw'
581            PRINT*,albedo(igout,1),emis(igout),mu0(igout),
582     s           fract(igout), fluxsurf_lw(igout),
583     $     fluxsurf_sw(igout,1)+fluxsurf_sw(igout,2)
584            PRINT*,'Tlay Tlev Play Plev dT/dt SW dT/dt LW (K/s)'
585            PRINT*,'daysec',daysec
586            DO l=1,nlayer
587               PRINT*,pt(igout,l),ptlev(igout,l),
588     s         pplay(igout,l),pplev(igout,l),
589     s         dtsw(igout,l),dtlw(igout,l)
590            ENDDO
591         ENDIF
592
593
594      END SUBROUTINE callradite
595
596      END MODULE callradite_mod
Note: See TracBrowser for help on using the repository browser.