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

Last change on this file since 1944 was 1918, checked in by emillour, 7 years ago

Mars GCM:
Code cleanup:

  • remove "comorbit.h" since it is no longer used.
  • turn "datafile.h" into module datafile_mod.F90 (and rename variable "datafile" as "datadir" since it stores the path to the datafile directory).

EM

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, only: aeropacity
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_mod. 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! compute ndomain
276! AS: moved out of firstcall to allow nesting+evoluting domain
277! ------------------------------------------------------------
278      ndomain= (ngrid-1) / ndomainsz + 1
279
280      IF (firstcall) THEN
281
282         write(*,*) 'Splitting radiative calculations: ',
283     $              ' ngrid,ndomainsz,ndomain',
284     $                ngrid,ndomainsz,ndomain
285
286c        Assign a number to the different scatterers
287c        -------------------------------------------
288
289         iaer_dust_conrath=0
290         iaer_dust_doubleq=0
291         iaer_dust_submicron=0
292         iaer_h2o_ice=0
293
294         aer_count=0
295         if (.NOT.active) then
296           do iaer=1,naerkind
297             if (name_iaer(iaer).eq."dust_conrath") then
298               iaer_dust_conrath = iaer
299               aer_count = aer_count + 1
300             endif
301           enddo
302         endif
303         if (doubleq.AND.active) then
304           do iaer=1,naerkind
305             if (name_iaer(iaer).eq."dust_doubleq") then
306               iaer_dust_doubleq = iaer
307               aer_count = aer_count + 1
308             endif
309           enddo
310         endif
311         if (submicron.AND.active) then
312           do iaer=1,naerkind
313             if (name_iaer(iaer).eq."dust_submicron") then
314               iaer_dust_submicron = iaer
315               aer_count = aer_count + 1
316             endif
317           enddo
318         endif
319         if (water.AND.activice) then
320           do iaer=1,naerkind
321             if (name_iaer(iaer).eq."h2o_ice") then
322               iaer_h2o_ice = iaer
323               aer_count = aer_count + 1
324             endif
325           enddo
326         endif
327
328c        Check that we identified all tracers:
329         if (aer_count.ne.naerkind) then
330           write(*,*) "callradite: found only ",aer_count," scatterers"
331           write(*,*) "               expected ",naerkind
332           write(*,*) "please make sure that the number of"
333           write(*,*) "scatterers in scatterers.h, the names"
334           write(*,*) "in callradite.F, and the flags in"
335           write(*,*) "callphys.def are all consistent!"
336           do iaer=1,naerkind
337             write(*,*)'      ',iaer,' ',trim(name_iaer(iaer))
338           enddo
339           stop
340         else
341           write(*,*) "callradite: found all scatterers, namely:"
342           do iaer=1,naerkind
343             write(*,*)'      ',iaer,' ',trim(name_iaer(iaer))
344           enddo
345         endif
346c        -------------------------------------------
347
348         gcp = g/cpp
349
350c        Loading the optical properties in external look-up tables:
351         CALL SUAER
352!         CALL SULW ! this step is now done in ini_yomlw_h
353
354         if (ngrid .EQ. 1) then
355           if (ndomainsz .NE. 1) then
356             print*
357             print*,'ATTENTION !!!'
358             print*,'pour tourner en 1D, '
359             print*,'fixer ndomainsz=1 dans phymars/dimradmars_h'
360             print*
361             call exit(1)
362           endif
363         endif
364
365         firstcall=.false.
366      END IF
367
368c     Computing aerosol optical properties and opacity
369c     ------------------------------------------------
370
371c     Updating aerosol size distributions:
372      CALL updatereffrad(ngrid,nlayer,
373     &                rdust,rice,nuice,
374     &                reffrad,nueffrad,
375     &                pq,tauscaling,tau,pplay)
376
377c     Computing 3D scattering parameters:
378      CALL aeroptproperties(ngrid,nlayer,reffrad,nueffrad,
379     &                      QVISsQREF3d,omegaVIS3d,gVIS3d,
380     &                      QIRsQREF3d,omegaIR3d,gIR3d,
381     &                      QREFvis3d,QREFir3d,
382     &                      omegaREFvis3d,omegaREFir3d)
383
384c     Computing aerosol optical depth in each layer:
385      CALL aeropacity(ngrid,nlayer,nq,zday,pplay,pplev,ls,
386     &     pq,tauscaling,tauref,tau,taucloudtes,aerosol,dsodust,reffrad,
387     &     nueffrad,QREFvis3d,QREFir3d,omegaREFvis3d,omegaREFir3d,
388     &     clearsky,totcloudfrac)
389
390c     Starting loop on sub-domain
391c     ----------------------------
392
393      DO jd=1,ndomain
394        ig0=(jd-1)*ndomainsz
395        if (jd.eq.ndomain) then
396         nd=ngrid-ig0
397        else
398         nd=ndomainsz
399        endif
400
401c       Spliting input variable in sub-domain input variables
402c       ---------------------------------------------------
403
404        do l=1,nlaylte
405         do ig = 1,nd
406           do iaer = 1, naerkind
407             do ich = 1, nsun
408               zQVISsQREF3d(ig,l,ich,iaer) =
409     &                           QVISsQREF3d(ig0+ig,l,ich,iaer)
410               zomegaVIS3d(ig,l,ich,iaer) =
411     &                           omegaVIS3d(ig0+ig,l,ich,iaer)
412               zgVIS3d(ig,l,ich,iaer) =
413     &                           gVIS3d(ig0+ig,l,ich,iaer)
414             enddo
415             do ich = 1, nir
416               zQIRsQREF3d(ig,l,ich,iaer) =
417     &                           QIRsQREF3d(ig0+ig,l,ich,iaer)
418               zomegaIR3d(ig,l,ich,iaer) =
419     &                           omegaIR3d(ig0+ig,l,ich,iaer)
420               zgIR3d(ig,l,ich,iaer) =
421     &                           gIR3d(ig0+ig,l,ich,iaer)
422             enddo
423           enddo
424         enddo
425        enddo
426
427        do l=1,nlaylte+1
428         do ig = 1,nd
429          zplev(ig,l) = pplev(ig0+ig,l)
430         enddo
431        enddo
432
433        do l=1,nlaylte
434         do ig = 1,nd
435          zplay(ig,l) = pplay(ig0+ig,l)
436          zt(ig,l) = pt(ig0+ig,l)
437c         Thickness of each layer (Pa) :
438          zdp(ig,l)= pplev(ig0+ig,l) - pplev(ig0+ig,l+1)
439         enddo
440        enddo
441
442        do n=1,naerkind
443          do l=1,nlaylte
444            do ig=1,nd
445              zaerosol(ig,l,n) = aerosol(ig0+ig,l,n)
446            enddo
447          enddo
448        enddo
449
450        do j=1,2
451          do ig = 1,nd
452           zalbedo(ig,j) = albedo(ig0+ig,j)
453          enddo
454        enddo
455
456c       Intermediate  levels: (computing tlev)
457c       ---------------------------------------
458c       Extrapolation for the air temperature above the surface
459        DO ig=1,nd
460              zztlev(ig,1)=zt(ig,1)+
461     s        (zplev(ig,1)-zplay(ig,1))*
462     s        (zt(ig,1)-zt(ig,2))/(zplay(ig,1)-zplay(ig,2))
463
464              zdt0(ig) = tsurf(ig0+ig) - zztlev(ig,1)
465        ENDDO
466
467        DO l=2,nlaylte
468         DO ig=1, nd
469               zztlev(ig,l)=0.5*(zt(ig,l-1)+zt(ig,l))
470         ENDDO
471        ENDDO
472
473        DO ig=1, nd
474           zztlev(ig,nlaylte+1)=zt(ig,nlaylte)
475        ENDDO
476
477
478c       Longwave ("lw") radiative transfer (= thermal infrared)
479c       -------------------------------------------------------
480        call lwmain (ig0,icount,nd,nflev
481     .        ,zdp,zdt0,emis(ig0+1),zplev,zztlev,zt
482     .        ,zaerosol,zzdtlw
483     .        ,fluxsurf_lw(ig0+1),fluxtop_lw(ig0+1)
484     .        ,znetrad
485     &        ,zQIRsQREF3d,zomegaIR3d,zgIR3d
486     &        ,co2ice(ig0+1))
487
488c       Shortwave ("sw") radiative transfer (= solar radiation)
489c       -------------------------------------------------------
490c          Mars solar constant (W m-2)
491c          1370 W.m-2 is the solar constant at 1 AU.
492           cste_mars=1370./(dist_sol*dist_sol)
493
494           call swmain ( nd, nflev,
495     S     cste_mars, zalbedo,
496     S     mu0(ig0+1), zdp, zplev, zaerosol, fract(ig0+1),
497     S     zzdtsw, zfluxd_sw, zfluxu_sw,
498     &     zQVISsQREF3d,zomegaVIS3d,zgVIS3d)
499
500c       ------------------------------------------------------------
501c       Un-spliting output variable from sub-domain input variables
502c       ------------------------------------------------------------
503
504        do l=1,nlaylte
505         do ig = 1,nd
506          dtlw(ig0+ig,l) = zzdtlw(ig,l)
507          dtsw(ig0+ig,l) = zzdtsw(ig,l)
508         enddo
509        enddo
510
511        do l=1,nlaylte+1
512         do ig = 1,nd
513          ptlev(ig0+ig,l) = zztlev(ig,l)
514         enddo
515        enddo
516
517        do ig = 1,nd
518          fluxsurf_sw(ig0+ig,1) = zfluxd_sw(ig,1,1)
519          fluxsurf_sw(ig0+ig,2) = zfluxd_sw(ig,1,2)
520          fluxtop_sw(ig0+ig,1) = zfluxu_sw(ig,nlaylte+1,1)
521          fluxtop_sw(ig0+ig,2) = zfluxu_sw(ig,nlaylte+1,2)
522        enddo
523
524      ENDDO         !   (boucle jd=1, ndomain)
525
526c     Zero tendencies for any remaining layers between nlaylte and nlayer
527      if (nlayer.gt.nlaylte) then
528         do l = nlaylte+1, nlayer
529            do ig = 1, ngrid
530               dtlw(ig, l) = 0.
531               dtsw(ig, l) = 0.
532            enddo
533         enddo
534      endif
535
536c     Output for debugging if lwrite=T
537c     --------------------------------
538c     Write all nlayer layers, even though only nlaylte layers may have
539c     non-zero tendencies.
540
541         IF(lwrite) THEN
542            PRINT*,'Diagnotique for the radiation'
543            PRINT*,'albedo, emissiv, mu0,fract,fluxsurf_lw,fluxsurf_sw'
544            PRINT*,albedo(igout,1),emis(igout),mu0(igout),
545     s           fract(igout), fluxsurf_lw(igout),
546     $     fluxsurf_sw(igout,1)+fluxsurf_sw(igout,2)
547            PRINT*,'Tlay Tlev Play Plev dT/dt SW dT/dt LW (K/s)'
548            PRINT*,'daysec',daysec
549            DO l=1,nlayer
550               PRINT*,pt(igout,l),ptlev(igout,l),
551     s         pplay(igout,l),pplev(igout,l),
552     s         dtsw(igout,l),dtlw(igout,l)
553            ENDDO
554         ENDIF
555
556
557      END SUBROUTINE callradite
558
559      END MODULE callradite_mod
Note: See TracBrowser for help on using the repository browser.