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

Last change on this file since 2154 was 1775, checked in by aslmd, 7 years ago

LMDZ.MARS setting the stage for maybe fixing nesting in the LMD_MM_MARS 4. ensure arrays not allocated in read_dust_scenario if already allocated -- other changes are useful comments for subsequent developments

File size: 1.8 KB
RevLine 
[38]1      subroutine nlthermeq(ngrid, nlayer, pplev, pplay)
2c
[1047]3c  Compute the number of layers nlaylte (stored in module yomlw_h)
[38]4c  over which local thermodynamic equilibrium
5c  radiation scheme should be run to be sure of covering at least to a
6c  height greater than (pressure lower than) p=pminte, set in nlteparams.h.
7c  The maximum layer needed is found for the worst possible case.
8c  Stephen Lewis 6/2000
9c  Modified Y. Wanherdrick/ F. Forget 09/2000
[1047]10      use yomlw_h, only: nlaylte
[38]11      implicit none
12#include "nlteparams.h"
13#include "callkeys.h"
14
15c
16c     Input:
17      integer ngrid, nlayer
18      real pplev(ngrid, nlayer+1)
19      real pplay(ngrid, nlayer)
20c
21c     Local:
22      integer igpmax, ismax
23      logical firstcall
24      data firstcall /.true./
25      save firstcall, igpmax
26c
27      if(firstcall) then
28c     Find the location of maximum surface pressure.
29c     Location won't vary much so only do it at the start;
30c     with no topography location would vary, but this is only
31c     needed for an estimate so any point would do in that case.
[1775]32!!    AS: can be problem w MESOSCALE nesting (ignored for the moment)
[38]33         igpmax = ismax(ngrid, pplev, 1)
34         write(*, 10) ptrans
35         write(*, 20) zw
36         write(*, 30) pminte
37         firstcall = .false.
38      endif
39c
40      IF(callnlte) then
41c       Find first layer above pminte at this location
42        do nlaylte = nlayer, 1, -1
43           if (pplay(igpmax, nlaylte).gt.pminte)  go to 100
44        enddo
45      ELSE
46        nlaylte=nlayer
47      END IF
48  100 write(*,*) 'LTE rad. calculations up to layer ',  nlaylte
49c
50      return
51c
52   10 format(' nlthermeq: transition to NLTE centred at ',f6.2,'Pa')
53   20 format('               half-width (scale heights) ',f6.2)
54   30 format('          suggested LTE coverage at least ',f6.2,'Pa')
55      end
Note: See TracBrowser for help on using the repository browser.