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

Last change on this file since 1655 was 1524, checked in by emillour, 9 years ago

All GCMS:
More updates to enforce dynamics/physics separation:

get rid of references to "temps_mod" from physics packages;
make a "time_phylmdz_mod.F90" module to store that
information and fill it via "iniphysiq".

EM

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