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

Last change on this file since 2347 was 2252, checked in by abierjon, 5 years ago

Mars GCM:
Bug fix following r2248 in aeropacity_mod and topmons_mod : since dsodust, dsords and dsotop are diagnostic physiq_mod variables, we don't want them to be reinitialized at each call of aeropacity_mod and topmons_mod, but we initialize them once and for all at the beginning of physiq_mod instead.
AB

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