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

Last change on this file since 1880 was 1774, checked in by aslmd, 7 years ago

LMDZ.MARS setting the stage for maybe fixing nesting in the LMD_MM_MARS 3. moving out domain-dependent computations from firstcalls. these are few lines with simple computations, no impact on runtime.

File size: 5.5 KB
Line 
1      SUBROUTINE calldrag_noro(ngrid,nlayer,ptimestep,
2     &                 pplay,pplev,pt,pu,pv,pdtgw,pdugw,pdvgw)
3
4
5
6       use surfdat_h, only: zstd, zsig, zgam, zthe
7       use dimradmars_mod, only: ndomainsz
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
19c        The sub-grid size is defined in dimradmars_mod.
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
74      REAL sigtest(nlayer+1)
75      INTEGER igwd,igwdim,itest(ngrid)
76
77      INTEGER :: ndomain
78!      parameter (ndomain = (ngrid-1) / ndomainsz + 1)
79
80      INTEGER l,ig
81      INTEGER jd,ig0,nd
82
83      REAL zulow(ngrid),zvlow(ngrid)
84      REAL zustr(ngrid),zvstr(ngrid)
85
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)
91      INTEGER zidx(ndomainsz)
92      REAL zzdtgw(ndomainsz,nlayer)
93      REAL zzdugw(ndomainsz,nlayer)
94      REAL zzdvgw(ndomainsz,nlayer)
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
113
114         do l=1,nlayer+1
115           sigtest(l)=pplev(1,l)/pplev(1,1)
116         enddo
117         call sugwd(nlayer,sigtest)
118
119         if (ngrid .EQ. 1) then
120           if (ndomainsz .NE. 1) then
121             print*
122             print*,'ATTENTION !!!'
123             print*,'pour tourner en 1D, meme pour drag_noro '
124             print*,'fixer ndomainsz=1 dans phymars/dimradmars_mod'
125             print*
126             call exit(1)
127           endif
128         endif
129
130         firstcall=.false.
131      END IF
132
133      !! AS: moved out of firstcall to allow nesting+evoluting horiz domain
134      ndomain = (ngrid-1) / ndomainsz + 1
135
136c     Starting loop on sub-domain
137c     ----------------------------
138
139      DO jd=1,ndomain
140        ig0=(jd-1)*ndomainsz
141        if (jd.eq.ndomain) then
142          nd=ngrid-ig0
143        else
144          nd=ndomainsz
145        endif
146
147c       Detecting points concerned by the scheme
148c       ----------------------------------------
149
150        igwd=0
151        DO ig=ig0+1,ig0+nd
152          itest(ig)=0
153          ll=zstd(ig).gt.50.0
154          IF(ll) then
155            itest(ig)=1
156            igwd=igwd+1
157            zidx(igwd)=ig - ig0
158          ENDIF
159        ENDDO
160        IGWDIM=MAX(1,IGWD)
161
162c       Spliting input variable in sub-domain input variables
163c       ---------------------------------------------------
164
165        do l=1,nlayer+1
166          do ig = 1,nd
167           zplev(ig,l) = pplev(ig0+ig,l)
168          enddo
169        enddo
170
171        do l=1,nlayer
172          do ig = 1,nd
173            zplay(ig,l) = pplay(ig0+ig,l)
174            zt(ig,l) = pt(ig0+ig,l)
175            zu(ig,l) = pu(ig0+ig,l)
176            zv(ig,l) = pv(ig0+ig,l)
177          enddo
178        enddo
179
180c       Calling gravity wave and subgrid scale topo parameterization
181c       -------------------------------------------------------------
182
183        call drag_noro (nd,nlayer,ptimestep,zplay,zplev,
184     e        zstd(ig0+1),zsig(ig0+1),zgam(ig0+1),zthe(ig0+1),
185     e        igwd,igwdim,zidx,itest(ig0+1),
186     e        zt, zu, zv,
187     s        zulow(ig0+1),zvlow(ig0+1),zustr(ig0+1),zvstr(ig0+1),
188     s        zzdtgw,zzdugw,zzdvgw)
189
190c       Un-spliting output variable from sub-domain input variables
191c       ------------------------------------------------------------
192c       (and devide by ptimestep -> true tendancies)
193
194        do l=1,nlayer
195         do ig = 1,nd
196          pdtgw(ig0+ig,l) = zzdtgw(ig,l)/ptimestep
197          pdugw(ig0+ig,l) = zzdugw(ig,l)/ptimestep
198          pdvgw(ig0+ig,l) = zzdvgw(ig,l)/ptimestep
199         enddo
200        enddo
201
202      ENDDO         !   (boucle jd=1, ndomain)
203
204      return
205      end
206
Note: See TracBrowser for help on using the repository browser.