source: trunk/LMDZ.MARS/libf/phymars/lwmain_mod.F @ 3026

Last change on this file since 3026 was 3004, checked in by emillour, 18 months ago

Mars PCM:
Some code cleanup. Remove obsolete "comg1d.h" and "writeg1d.F" (were used to
specifically output for Grads in 1D).
Also turned lwi and lwflux into modules while at it.
EM

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