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