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

Last change on this file since 2616 was 2584, checked in by romain.vande, 3 years ago

Second stage of implementation of Open_MP in the physic.
Run with callrad=.true.

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