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

Last change on this file since 4063 was 3901, checked in by emillour, 6 months ago

Mars PCM:
Some code tidying:

  • turn aeroptproperties, albedocaps, cvmgp and convadj into modules
  • remove useless check in callradite
  • clean module vlz_fi (remove obsolete #ifdef CRAY cpp alternatives)
  • remove function cvmgp since it was only called under #ifdef CRAY in vlz_fi

EM

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