source: trunk/LMDZ.VENUS/libf/phyvenus/nlthermeq.F @ 3094

Last change on this file since 3094 was 1591, checked in by slebonnois, 8 years ago

SL: Mise a jour de la haute atmosphere, du transfert radiatif (solaire=Rainer Haus; IR: reglages de juin 2016), et implementation des variations de temperature potentielle dans la basse atmosphere (variation de la masse molaire moyenne)

File size: 2.0 KB
Line 
1      subroutine nlthermeq(nlon, nlayer, pplev, pplay)
2c
3c  Compute the number of layers nlaylte (stored in common 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 dimphy
11      implicit none
12c#include "dimradmars.h"
13#include "nlteparams.h"
14c#include "yomlw.h"
15#include "clesphys.h"
16
17c
18c     Input:
19      integer nlon, nlayer
20      real pplev(nlon, nlayer+1)
21      real pplay(nlon, nlayer)
22c
23c     Local:
24      integer igpmax, ismax
25      logical firstcall
26      data firstcall /.true./
27      save firstcall, igpmax
28
29
30      INTEGER i,ix
31      real sxmax
32
33ccc     
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      ismax=1 
40      sxmax=pplev(1,1)
41      ix=1
42        do i=1,nlon-1
43         if(pplev(i,ix).gt.sxmax) then
44           sxmax=pplev(i,ix)
45           ismax=i+1
46         endif
47       enddo
48
49
50         igpmax = ismax            ! longitude/ latitude where pression is maximum
51         write(*, 10) ptrans
52         write(*, 20) zw
53         write(*, 30) pminte
54         firstcall = .false.
55      endif
56
57      IF(callnlte .or. callnirco2) THEN
58c       Find first layer above pminte at this location
59        do nlaylte = nlayer, 1, -1
60      if (pplay(igpmax, nlaylte).gt.pminte)  go to 100
61        enddo
62      ELSE
63        nlaylte=nlayer       
64      END IF
65
66c
67 100    return
68c
69
70   10 format(' nlthermeq: transition to NLTE centred at ',f6.2,'Pa')
71   20 format('               half-width (scale heights) ',f6.2)
72   30 format('          suggested LTE coverage at least ',f6.2,'Pa')
73   40 format(' nlthermeq: purely NLTE contribution over (nlayer) ',f6.4)
74
75      end
Note: See TracBrowser for help on using the repository browser.