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

Last change on this file since 3807 was 3726, checked in by emillour, 2 months ago

Mars PCM:
Turn "callkeys.h" into module "callkeys_mod.F90"
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      use callkeys_mod, only: callnlte
19      implicit none
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 ismax
29      integer,save :: igpmax
30      logical,save :: firstcall = .true.
31!$OMP THREADPRIVATE(firstcall,igpmax)
32
33c
34      if(firstcall) then
35c     Find the location of maximum surface pressure.
36c     Location won't vary much so only do it at the start;
37c     with no topography location would vary, but this is only
38c     needed for an estimate so any point would do in that case.
39!!    AS: can be problem w MESOSCALE nesting (ignored for the moment)
40!! Ehouarn: this could also be a problem when in parallel with the GCM as well...
41!! fortunately with current setup (pminte=0.04) and hybrid vertical coordinates
42!! the transition point is in purely pressure layers and thus always the same.
43!! To be cleaned and made more general and robust someday...
44         igpmax = ismax(ngrid, pplev, 1)
45         write(*, 10) ptrans
46         write(*, 20) zw
47         write(*, 30) pminte
48         firstcall = .false.
49      endif
50c
51      IF(callnlte) then
52c       Find first layer above pminte at this location
53        do nlaylte = nlayer, 1, -1
54           if (pplay(igpmax, nlaylte).gt.pminte)  go to 100
55        enddo
56      ELSE
57        nlaylte=nlayer
58      END IF
59  100 write(*,*) 'LTE rad. calculations up to layer ',  nlaylte
60c
61      return
62c
63   10 format(' nlthermeq: transition to NLTE centred at ',f6.2,'Pa')
64   20 format('               half-width (scale heights) ',f6.2)
65   30 format('          suggested LTE coverage at least ',f6.2,'Pa')
66      end subroutine nlthermeq
67
68      END MODULE nlthermeq_mod
Note: See TracBrowser for help on using the repository browser.