source: trunk/LMDZ.MARS/libf/phymars/lwmain.F @ 1233

Last change on this file since 1233 was 1047, checked in by emillour, 11 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File size: 7.5 KB
Line 
1       subroutine lwmain (ig0,icount,kdlon,kflev
2     .                   ,dp,dt0,emis
3     .                   ,plev,tlev,tlay,aerosol,coolrate
4     .                   ,fluxground,fluxtop
5     .                   ,netrad
6     &                   ,QIRsQREF3d,omegaIR3d,gIR3d
7     &                   ,co2ice)
8
9c----------------------------------------------------------------------
10c     LWMAIN     organizes the LTE longwave calculations
11c     for layer 1 to layer "nlaylte" (stored in  "yomlw_h")
12c----------------------------------------------------------------------
13
14      use dimradmars_mod, only: ndlo2, nflev, nir, ndlon, nuco2
15      use yomlw_h, only: nlaylte, xi
16      implicit none
17 
18!#include "dimensions.h"
19!#include "dimphys.h"
20!#include "dimradmars.h"
21#include "callkeys.h"
22#include "comg1d.h"
23! naerkind is set in scatterers.h (built when compiling with makegcm -s #)
24#include"scatterers.h"
25!#include "yomlw.h"
26
27c----------------------------------------------------------------------
28c         0.1   arguments
29c               ---------
30c                                                            inputs:
31c                                                            -------
32      integer ig0
33      integer icount
34      integer kdlon            ! part of ngrid
35      integer kflev            ! part of nlayer
36
37      real dp (ndlo2,kflev)         ! layer pressure thickness (Pa)
38      real dt0 (ndlo2)              ! surface temperature discontinuity (K)
39      real emis (ndlo2)             ! surface emissivity
40      real plev (ndlo2,kflev+1)     ! level pressure (Pa)
41      real tlev (ndlo2,kflev+1)     ! level temperature (K)
42      real tlay (ndlo2,kflev)       ! layer temperature (K)
43      real aerosol(ndlo2,kflev,naerkind)      !  aerosol extinction optical
44c                         depth at reference wavelength "longrefvis" set
45c                         in dimradmars_mod , in each layer, for one of
46c                         the "naerkind" kind of aerosol optical properties.
47
48
49c                                                            outputs:
50c                                                            --------
51      real coolrate(ndlo2,kflev)      ! cooling rate (K/s)
52      real fluxground(ndlo2)          ! downward ground flux (W/m2)
53      real fluxtop(ndlo2)             ! outgoing upward flux (W/m2) ("OLR")
54      real netrad (ndlo2,kflev)       ! radiative budget (W/m2)
55c     Aerosol optical properties
56      REAL :: QIRsQREF3d(ndlo2,kflev,nir,naerkind)
57      REAL :: omegaIR3d(ndlo2,kflev,nir,naerkind)
58      REAL :: gIR3d(ndlo2,kflev,nir,naerkind)
59
60c----------------------------------------------------------------------
61c         0.2   local arrays
62c               ------------
63
64      real aer_t (ndlon,nuco2,nflev+1)  ! transmission (aer)
65      real co2_u (ndlon,nuco2,nflev+1)  ! absorber amounts (co2)
66      real co2_up (ndlon,nuco2,nflev+1) ! idem scaled by the pressure (co2)
67
68      real bsurf (ndlon,nir)            ! surface spectral planck function
69      real btop (ndlon,nir)             ! top spectral planck function
70      real blev (ndlon,nir,nflev+1)     ! level   spectral planck function
71      real blay (ndlon,nir,nflev)       ! layer   spectral planck function
72      real dblay (ndlon,nir,nflev)      ! layer gradient spectral planck function
73      real dbsublay (ndlon,nir,2*nflev) ! layer gradient spectral planck function
74                                        ! in sub layers
75
76      real tautotal(ndlon,nflev,nir)  ! \   Total single scattering
77      real omegtotal(ndlon,nflev,nir) !  >  properties (Addition of the
78      real gtotal(ndlon,nflev,nir)    ! /   NAERKIND aerosols prop.)
79
80      real newcoolrate(ndlon,nflev) ! cooling rate (K/s) / with implicite scheme
81
82      REAL co2ice(ndlo2)           ! co2 ice surface layer (kg.m-2)
83      REAL emis_gaz(ndlo2)         ! emissivity for gaz computations
84
85      integer jk,jkk,ja,jl
86
87      logical firstcall
88      save firstcall
89      data firstcall/.true./
90
91
92c----------------------------------------------------------------------
93c         0.3   Initialisation
94c               --------------
95
96      if (firstcall) then
97
98        firstcall = .false.
99
100        xi (:,:,:,:)=0.
101
102      endif
103
104      DO jl=1 , kdlon
105         IF(co2ice(jl) .GT. 20.e-3) THEN
106             emis_gaz(jl)=1.
107         ELSE
108             emis_gaz(jl)=emis(jl)
109         ENDIF
110      ENDDO
111
112c----------------------------------------------------------------------
113c         1.0   planck function
114c               ---------------
115
116      call lwb ( kdlon, kflev, tlev, tlay, dt0
117     .         , bsurf, btop, blay, blev, dblay, dbsublay)
118
119c----------------------------------------------------------------------
120c         2.0   absorber amounts
121c               ----------------
122
123      call lwu ( kdlon, kflev
124     .         , dp, plev, tlay, aerosol
125     &         , QIRsQREF3d,omegaIR3d,gIR3d
126     .         , aer_t, co2_u, co2_up
127     .         , tautotal,omegtotal,gtotal)
128
129c----------------------------------------------------------------------
130c         3.0   transmission functions / exchange coefficiants
131c               ----------------------------------------------
132
133c                                                                distants
134c                                                                --------
135                    if( mod(icount-1,ilwd).eq.0) then
136
137c     print*, 'CALL of DISTANTS'
138      call lwxd ( ig0, kdlon, kflev, emis_gaz
139     .          , aer_t, co2_u, co2_up)
140
141                    endif
142c                                                              neighbours
143c                                                              ----------
144                    if( mod(icount-1,ilwn).eq.0) then
145
146c     print*, 'CALL of NEIGHBOURS'
147      call lwxn ( ig0, kdlon, kflev
148     .          , dp
149     .          , aer_t, co2_u, co2_up)
150
151                    endif
152c                                                              boundaries
153c                                                              ----------
154                    if( mod(icount-1,ilwb).eq.0) then
155
156c     print*, 'CALL of BOUNDARIES'
157      call lwxb ( ig0, kdlon, kflev, emis_gaz
158     .          , aer_t, co2_u, co2_up)
159
160                    endif
161
162c----------------------------------------------------------------------
163c         4.0   cooling rate
164c               ------------
165
166      call lwflux ( ig0, kdlon, kflev, dp
167     .            , bsurf, btop, blev, blay, dbsublay
168     .            , tlay, tlev, dt0      ! pour sortie dans g2d uniquement
169     .            , emis
170     .            , tautotal,omegtotal,gtotal
171     .            , coolrate, fluxground, fluxtop
172     .            , netrad)
173
174c     do jk = 1, nlaylte
175c       print*,coolrate(1,jk)
176c     enddo
177     
178c       do jkk = 0 , nlaylte+1
179c         do jk = 0 , nlaylte+1
180c           do ja = 1 , nuco2
181c             do jl = 1 , ngrid
182c      if (xi (jl,ja,jk,jkk) .LT. 0
183c    .       .OR. xi (jl,ja,jk,jkk) .GT. 1 ) then
184c                 print*,'xi bande',ja,jk,jkk,xi (jl,ja,jk,jkk)
185c      endif
186c             enddo
187c           enddo
188c         enddo
189c       enddo
190
191c----------------------------------------------------------------------
192c
193c          5.    shema semi-implicite  (lwi)
194c                ---------------------------
195c
196c
197      call lwi (ig0,kdlon,kflev,netrad,dblay,dp
198     .          , newcoolrate)
199c
200c  Verif que   (X sol,space) + somme(X i,sol) = 1
201c
202      do jkk = 1 , nlaylte
203        do jl = 1 , kdlon
204c     print*,'NEW et OLD coolrate :',jkk,newcoolrate(jl,jkk)
205c    .  ,coolrate(jl,jkk)
206      coolrate(jl,jkk) = newcoolrate(jl,jkk)
207        enddo
208      enddo
209c
210c----------------------------------------------------------------------
211
212      return
213      end
Note: See TracBrowser for help on using the repository browser.