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

Last change on this file since 259 was 172, checked in by emillour, 13 years ago

Mars GCM

  • minor correction to callradite.F (to enable compilation in debug mode with ifort when there is only one tracer).

EM

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