source: trunk/LMDZ.MARS/libf/phymars/callradite.F @ 363

Last change on this file since 363 was 358, checked in by aslmd, 14 years ago

LMDZ.MARS:


A FIRST VERSION WITH SCAVENGING
TRANSPARENT TO CASUAL USER
FOR WATER CYCLE FOLKS, SEE BELOW


[TN and AS solved conflicts due to r330 >> r357]


07/11/11 == JBM

Changed watercloud.F to call two separate routines,

simpleclouds.F or improvedclouds.F, which are a simplified and
full microphysical scheme for cloud formation, respectively.
Removed the tag called "improved" in watercloud.F, and added
another tag called "microphys" which is defined in callphys.F
instead. Changed routines: callkeys, inifis, physiq, watercloud.

Reimplemented the use of typical profiles of dust particle sizes

and CCNs in simpleclouds.F. Corrected the previously used
analytical CCN profile. Moved ccn_factor to simpleclouds.F,
which won't be used in the new microphysical scheme. Changed
routines: aeropacity, physiq, simpleclouds, watercloud.

Computed rdust at the end of callsedim instead of updatereffrad,

except at firstcall. Renamed rfall into rsedcloud and moved it
in simpleclouds. Moved nuice_sed to initracer.F and added it to
tracer.h. Changed routines: callsedim, physiq, tracer.h,
watercloud, initracer, simpleclouds, updatereffrad.

Added two tracers for the CCNs. Added the different tests in

initracer.F to make sure that, for example, microphys is not
used without doubleq, etc. Corrected an inconsistency in
callsedim.F, and changed the way r0 is updated. Changes
routines: callsedim, inifis, initracer, physiq, testphys1d,
tracer.h.

Added the ability to have a spatially variable density in

newsedim (same method as that used for the radius of
sedimentation). Required the addition of one input to newsedim,
which is the size of the density variable "rho". Changed
routines: callsedim, newsedim.

Added an output to aeropacity called "tauscaling", which is a

factor to convert dust_mass and dust_number to true quantities,
based on tauvis. Removed ccn and qdust from aeropacity, which
are obsolete now.

Wrote improvedcloud.F which includes all the microphysical

scheme, and connected it to the sedimentation scheme. Added and
changed routines: callsedim, physiq, growthrate, nuclea,
improvedclouds, initracer, watercloud, microphys.h.

07/11/11 == TN

Added CCN & dust tracers updates in physiq.F

Corrected a bug that can give negative CCN numbers, due to the
use of single precision values (only 6 decimals) whereas up to 10E+10
CCN can be formed or disappears...

Corrected physical bug that causes h2o_vap or h2o_ice

to be negative in improvedclouds.F.

Corrected physical bug that causes CCN & dust tracers

to be negative in improvedclouds.F when all ice sublimates and
releases dust

Added parameter contact mteta in callphys.def

Default value is still 0.95, see inifis.F

Changed tendancy computation for dust_number in improvedclouds.F

that was not the right one. Indeed, the scavenged dust_number tracer
is computed from the dust_mass one, and its tendancy before scavenging
must be taken into account to compute its scavenging's tendancy.

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