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

Last change on this file since 3026 was 3016, checked in by emillour, 17 months ago

Mars PCM:
Code cleanup: nltedata.h is only included in nltecool.F so turn this data
into module data there. Also convert lteparams.h into module nlteparams_h.F90.
And while at it also turn nlthermeq.F and blendrad.F into modules.
EM

File size: 2.3 KB
Line 
1      MODULE nlthermeq_mod
2     
3      IMPLICIT NONE
4     
5      CONTAINS
6
7      subroutine nlthermeq(ngrid, nlayer, pplev, pplay)
8c
9c  Compute the number of layers nlaylte (stored in module yomlw_h)
10c  over which local thermodynamic equilibrium
11c  radiation scheme should be run to be sure of covering at least to a
12c  height greater than (pressure lower than) p=pminte, set in nlteparams.h.
13c  The maximum layer needed is found for the worst possible case.
14c  Stephen Lewis 6/2000
15c  Modified Y. Wanherdrick/ F. Forget 09/2000
16      use yomlw_h, only: nlaylte
17      use nlteparams_h, only: ptrans, pminte, zw
18      implicit none
19      include "callkeys.h"
20
21c
22c     Input:
23      integer,intent(in) :: ngrid, nlayer
24      real,intent(in) :: pplev(ngrid, nlayer+1)
25      real,intent(in) :: pplay(ngrid, nlayer)
26c
27c     Local:
28      integer igpmax, ismax
29      logical firstcall
30
31!$OMP THREADPRIVATE(firstcall,igpmax)
32
33      data firstcall /.true./
34      save firstcall, igpmax
35c
36      if(firstcall) then
37c     Find the location of maximum surface pressure.
38c     Location won't vary much so only do it at the start;
39c     with no topography location would vary, but this is only
40c     needed for an estimate so any point would do in that case.
41!!    AS: can be problem w MESOSCALE nesting (ignored for the moment)
42!! Ehouarn: this could also be a problem when in parallel with the GCM as well...
43!! fortunately with current setup (pminte=0.04) and hybrid vertical coordinates
44!! the transition point is in purely pressure layers and thus always the same.
45!! To be cleaned and made more general and robust someday...
46         igpmax = ismax(ngrid, pplev, 1)
47         write(*, 10) ptrans
48         write(*, 20) zw
49         write(*, 30) pminte
50         firstcall = .false.
51      endif
52c
53      IF(callnlte) then
54c       Find first layer above pminte at this location
55        do nlaylte = nlayer, 1, -1
56           if (pplay(igpmax, nlaylte).gt.pminte)  go to 100
57        enddo
58      ELSE
59        nlaylte=nlayer
60      END IF
61  100 write(*,*) 'LTE rad. calculations up to layer ',  nlaylte
62c
63      return
64c
65   10 format(' nlthermeq: transition to NLTE centred at ',f6.2,'Pa')
66   20 format('               half-width (scale heights) ',f6.2)
67   30 format('          suggested LTE coverage at least ',f6.2,'Pa')
68      end subroutine nlthermeq
69
70      END MODULE nlthermeq_mod
Note: See TracBrowser for help on using the repository browser.