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

Last change on this file since 1775 was 1774, checked in by aslmd, 8 years ago

LMDZ.MARS setting the stage for maybe fixing nesting in the LMD_MM_MARS 3. moving out domain-dependent computations from firstcalls. these are few lines with simple computations, no impact on runtime.

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