source: trunk/LMDZ.MARS/libf/phymars/calldrag_noro_mod.F @ 2090

Last change on this file since 2090 was 1912, checked in by emillour, 7 years ago

Mars GCM:
Tidying the gravity wave routines by turning them into modules:
orodrag.F -> orodrag_mod.F : note that the declared size of pvar(), which is
used in call to gwstress was wrong.
calldrag_noro.F -> calldrag_noro_mod.F
drag_noro.F -> drag_noro_mod.F
gwstress.F -> gwstress_mod.F
EM

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