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

Last change on this file since 1242 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: 1.8 KB
Line 
1      subroutine nlthermeq(ngrid, nlayer, pplev, pplay)
2c
3c  Compute the number of layers nlaylte (stored in module yomlw_h)
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
10      use yomlw_h, only: nlaylte
11      implicit none
12!#include "dimensions.h"
13!#include "dimphys.h"
14!#include "dimradmars.h"
15#include "nlteparams.h"
16!#include "yomlw.h"
17#include "callkeys.h"
18
19c
20c     Input:
21      integer ngrid, nlayer
22      real pplev(ngrid, nlayer+1)
23      real pplay(ngrid, nlayer)
24c
25c     Local:
26      integer igpmax, ismax
27      logical firstcall
28      data firstcall /.true./
29      save firstcall, igpmax
30c
31      if(firstcall) then
32c     Find the location of maximum surface pressure.
33c     Location won't vary much so only do it at the start;
34c     with no topography location would vary, but this is only
35c     needed for an estimate so any point would do in that case.
36         igpmax = ismax(ngrid, pplev, 1)
37         write(*, 10) ptrans
38         write(*, 20) zw
39         write(*, 30) pminte
40         firstcall = .false.
41      endif
42c
43      IF(callnlte) then
44c       Find first layer above pminte at this location
45        do nlaylte = nlayer, 1, -1
46           if (pplay(igpmax, nlaylte).gt.pminte)  go to 100
47        enddo
48      ELSE
49        nlaylte=nlayer
50      END IF
51  100 write(*,*) 'LTE rad. calculations up to layer ',  nlaylte
52c
53      return
54c
55   10 format(' nlthermeq: transition to NLTE centred at ',f6.2,'Pa')
56   20 format('               half-width (scale heights) ',f6.2)
57   30 format('          suggested LTE coverage at least ',f6.2,'Pa')
58      end
Note: See TracBrowser for help on using the repository browser.