source: trunk/LMDZ.MARS/libf/phymars/calldrag_noro.F @ 1711

Last change on this file since 1711 was 1266, checked in by aslmd, 11 years ago

LMDZ.MARS
IMPORTANT CHANGE

  • Remove all reference/use of nlayermx and dimphys.h
  • Made use of automatic arrays whenever arrays are needed with dimension nlayer
  • Remove lots of obsolete reference to dimensions.h
  • Converted iono.h and param_v4.h into corresponding modules

(with embedded subroutine to allocate arrays)
(no arrays allocated if thermosphere not used)

  • Deleted param.h and put contents into module param_v4_h
  • Adapted testphys1d, newstart, etc...
  • Made DATA arrays in param_read to be initialized by subroutine

fill_data_thermos in module param_v4_h

  • Optimized computations in paramfoto_compact (twice less dlog10 calculations)
  • Checked consistency before/after modification in debug mode
  • Checked performance is not impacted (same as before)
File size: 5.5 KB
RevLine 
[38]1      SUBROUTINE calldrag_noro(ngrid,nlayer,ptimestep,
2     &                 pplay,pplev,pt,pu,pv,pdtgw,pdugw,pdvgw)
3
4
5
[1047]6       use surfdat_h, only: zstd, zsig, zgam, zthe
7       use dimradmars_mod, only: ndomainsz
[38]8       IMPLICIT NONE
9c=======================================================================
10c   subject:
11c   --------
12c   Subroutine designed to call SUBROUTINE drag_noro
13c   Interface for sub-grid scale orographic scheme
14c   The purpose of this subroutine is
15c      1) Make some initial calculation at first call
16c      2) Split the calculation in several sub-grid
17c        ("sub-domain") to save memory and
18c        be able run on a workstation at high resolution
[1047]19c        The sub-grid size is defined in dimradmars_mod.
[38]20c
21c   author:   
22c   ------
23c           Christophe Hourdin/ Francois Forget
24c
25c   changes:
26c   -------
27c   > J.-B. Madeleine 10W12
28c   This version uses the variable's splitting, which can be usefull
29c     when performing very high resolution simulation like LES.
30c
31c   input:
32c   -----
33c   ngrid                 number of gridpoint of horizontal grid
34c   nlayer                Number of layer
35c   ptimestep             Physical timestep (s)
36c   pplay(ngrid,nlayer)    pressure (Pa) in the middle of each layer
37c   pplev(ngrid,nlayer+1)  pressure (Pa) at boundaries of each layer
38c   pt(ngrid,nlayer)       atmospheric temperature  (K)
39c   pu(ngrid,nlayer)       zonal wind (m s-1)
40c   pv(ngrid,nlayer)       meridional wind (m s-1)
41c
42c   output:
43c   -------
44c   pdtgw(ngrid,nlayer)    Temperature trend (K.s-1)
45c   pdugw(ngrid,nlayer)    zonal wind trend  (m.s-2)
46c   pdvgw(ngrid,nlayer)    meridional wind trend  (m.s-2)
47c
48c
49c
50c
51c
52c=======================================================================
53c
54c    0.  Declarations :
55c    ------------------
56c
57
58c-----------------------------------------------------------------------
59c    Input/Output
60c    ------------
61      INTEGER ngrid,nlayer 
62
63      real ptimestep
64
65      REAL pplev(ngrid,nlayer+1),pplay(ngrid,nlayer)
66      REAL pt(ngrid,nlayer), pu(ngrid,nlayer),pv(ngrid,nlayer)
67      REAL pdtgw(ngrid,nlayer), pdugw(ngrid,nlayer),pdvgw(ngrid,nlayer)
68
69
70c
71c    Local variables :
72c    -----------------
73
[1047]74      REAL sigtest(nlayer+1)
75      INTEGER igwd,igwdim,itest(ngrid)
[38]76
[1047]77      INTEGER,SAVE :: ndomain
78!      parameter (ndomain = (ngrid-1) / ndomainsz + 1)
[38]79
80      INTEGER l,ig
81      INTEGER jd,ig0,nd
82
[1047]83      REAL zulow(ngrid),zvlow(ngrid)
84      REAL zustr(ngrid),zvstr(ngrid)
[38]85
[1047]86      REAL zplev(ndomainsz,nlayer+1)
87      REAL zplay(ndomainsz,nlayer)
88      REAL zt(ndomainsz,nlayer)
89      REAL zu(ndomainsz,nlayer)
90      REAL zv(ndomainsz,nlayer)
[38]91      INTEGER zidx(ndomainsz)
[1047]92      REAL zzdtgw(ndomainsz,nlayer)
93      REAL zzdugw(ndomainsz,nlayer)
94      REAL zzdvgw(ndomainsz,nlayer)
[38]95
96      logical ll
97
98
99c   local saved variables
100c   ---------------------
101
102      LOGICAL firstcall
103      DATA firstcall/.true./
104      SAVE firstcall
105
106
107c----------------------------------------------------------------------
108
109c     Initialisation
110c     --------------
111
112      IF (firstcall) THEN
[1047]113        ndomain = (ngrid-1) / ndomainsz + 1
114
115         do l=1,nlayer+1
[38]116           sigtest(l)=pplev(1,l)/pplev(1,1)
117         enddo
[1047]118         call sugwd(nlayer,sigtest)
[38]119
[1047]120         if (ngrid .EQ. 1) then
[38]121           if (ndomainsz .NE. 1) then
122             print*
123             print*,'ATTENTION !!!'
124             print*,'pour tourner en 1D, meme pour drag_noro '
[1047]125             print*,'fixer ndomainsz=1 dans phymars/dimradmars_mod'
[38]126             print*
127             call exit(1)
128           endif
129         endif
130
131         firstcall=.false.
132      END IF
133
134c     Starting loop on sub-domain
135c     ----------------------------
136
137      DO jd=1,ndomain
138        ig0=(jd-1)*ndomainsz
139        if (jd.eq.ndomain) then
[1047]140          nd=ngrid-ig0
[38]141        else
142          nd=ndomainsz
143        endif
144
145c       Detecting points concerned by the scheme
146c       ----------------------------------------
147
148        igwd=0
149        DO ig=ig0+1,ig0+nd
150          itest(ig)=0
151          ll=zstd(ig).gt.50.0
152          IF(ll) then
153            itest(ig)=1
154            igwd=igwd+1
155            zidx(igwd)=ig - ig0
156          ENDIF
157        ENDDO
158        IGWDIM=MAX(1,IGWD)
159
160c       Spliting input variable in sub-domain input variables
161c       ---------------------------------------------------
162
163        do l=1,nlayer+1
164          do ig = 1,nd
165           zplev(ig,l) = pplev(ig0+ig,l)
166          enddo
167        enddo
168
169        do l=1,nlayer
170          do ig = 1,nd
171            zplay(ig,l) = pplay(ig0+ig,l)
172            zt(ig,l) = pt(ig0+ig,l)
173            zu(ig,l) = pu(ig0+ig,l)
174            zv(ig,l) = pv(ig0+ig,l)
175          enddo
176        enddo
177
178c       Calling gravity wave and subgrid scale topo parameterization
179c       -------------------------------------------------------------
180
181        call drag_noro (nd,nlayer,ptimestep,zplay,zplev,
182     e        zstd(ig0+1),zsig(ig0+1),zgam(ig0+1),zthe(ig0+1),
183     e        igwd,igwdim,zidx,itest(ig0+1),
184     e        zt, zu, zv,
185     s        zulow(ig0+1),zvlow(ig0+1),zustr(ig0+1),zvstr(ig0+1),
186     s        zzdtgw,zzdugw,zzdvgw)
187
188c       Un-spliting output variable from sub-domain input variables
189c       ------------------------------------------------------------
190c       (and devide by ptimestep -> true tendancies)
191
192        do l=1,nlayer
193         do ig = 1,nd
194          pdtgw(ig0+ig,l) = zzdtgw(ig,l)/ptimestep
195          pdugw(ig0+ig,l) = zzdugw(ig,l)/ptimestep
196          pdvgw(ig0+ig,l) = zzdvgw(ig,l)/ptimestep
197         enddo
198        enddo
199
200      ENDDO         !   (boucle jd=1, ndomain)
201
202      return
203      end
204
Note: See TracBrowser for help on using the repository browser.