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

Last change on this file since 1448 was 1353, checked in by tnavarro, 10 years ago

New option dustiropacity in callphys.def to change the reference IR opacity of dust + New output dsodust (density-scaled opacity). Without the use of this option, nothing changes for the uninformed user.

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