1 | !*********************************************************************** |
---|
2 | !* Copyright 2012,2013 * |
---|
3 | !* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, * |
---|
4 | !* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,* |
---|
5 | !* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, * |
---|
6 | !* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, * |
---|
7 | !* * |
---|
8 | !* This file is part of FLEXPART WRF * |
---|
9 | !* * |
---|
10 | !* FLEXPART is free software: you can redistribute it and/or modify * |
---|
11 | !* it under the terms of the GNU General Public License as published by* |
---|
12 | !* the Free Software Foundation, either version 3 of the License, or * |
---|
13 | !* (at your option) any later version. * |
---|
14 | !* * |
---|
15 | !* FLEXPART is distributed in the hope that it will be useful, * |
---|
16 | !* but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
17 | !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
18 | !* GNU General Public License for more details. * |
---|
19 | !* * |
---|
20 | !* You should have received a copy of the GNU General Public License * |
---|
21 | !* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
22 | !*********************************************************************** |
---|
23 | |
---|
24 | subroutine outgrid_init_reg() |
---|
25 | !******************************************************************************* |
---|
26 | ! * |
---|
27 | ! Note: This is the FLEXPART_WRF version of subroutine outgrid_init. * |
---|
28 | ! The computational grid is the WRF x-y grid rather than lat-lon. * |
---|
29 | ! * |
---|
30 | ! This routine calculates, for each grid cell of the output grid, the * |
---|
31 | ! volume, the surface area, and the areas of the northward and eastward * |
---|
32 | ! facing surfaces. * |
---|
33 | ! * |
---|
34 | ! Author: A. Stohl * |
---|
35 | ! * |
---|
36 | ! 7 August 2002 * |
---|
37 | ! * |
---|
38 | ! 26 Oct 2005, R. Easter - changes in gridarea, areaeast, areanorth * |
---|
39 | ! associated with WRF horizontal grid. * |
---|
40 | ! Dec 2005, R. Easter - changed names of "*lon0*" & "*lat0*" variables * |
---|
41 | ! July 2012, J. Brioude - modified for regular output grid * |
---|
42 | ! test version for https://www.flexpart.eu/ticket/307 |
---|
43 | !******************************************************************************* |
---|
44 | ! * |
---|
45 | ! Variables: * |
---|
46 | ! * |
---|
47 | ! area surface area of all output grid cells * |
---|
48 | ! areaeast eastward facing wall area of all output grid cells * |
---|
49 | ! areanorth northward facing wall area of all output grid cells * |
---|
50 | ! volume volumes of all output grid cells * |
---|
51 | ! * |
---|
52 | !******************************************************************************* |
---|
53 | |
---|
54 | use flux_mod |
---|
55 | use oh_mod |
---|
56 | use unc_mod |
---|
57 | use outg_mod |
---|
58 | use par_mod |
---|
59 | use com_mod |
---|
60 | |
---|
61 | ! include 'includepar' |
---|
62 | ! include 'includecom' |
---|
63 | implicit none |
---|
64 | integer :: ix,jy,kz,k,i,nage,l,iix,jjy,ixp,jyp,i1,j1,j,ngrid |
---|
65 | real :: ylat,gridarea,ylatp,ylatm,hzone,cosfact,cosfactm,cosfactp |
---|
66 | real :: ymet,xlon |
---|
67 | real :: xmet,xl,yl,ddx,ddy,rddx,rddy,p1,p2,p3,p4,xtn,ytn,oroh |
---|
68 | integer :: ks,kp,stat |
---|
69 | real,parameter :: eps=nxmax/3.e5 |
---|
70 | |
---|
71 | |
---|
72 | |
---|
73 | ! Compute surface area and volume of each grid cell: area, volume; |
---|
74 | ! and the areas of the northward and eastward facing walls: areaeast, areanorth |
---|
75 | !*********************************************************************** |
---|
76 | |
---|
77 | do jy=0,numygrid-1 |
---|
78 | |
---|
79 | ylat=outlat0+(float(jy)+0.5)*dyoutl |
---|
80 | ylatp=ylat+0.5*dyoutl |
---|
81 | ylatm=ylat-0.5*dyoutl |
---|
82 | if ((ylatm.lt.0).and.(ylatp.gt.0.)) then |
---|
83 | hzone=dyoutl*r_earth*pi180 |
---|
84 | else |
---|
85 | cosfact=cos(ylat*pi180)*r_earth |
---|
86 | cosfactp=cos(ylatp*pi180)*r_earth |
---|
87 | cosfactm=cos(ylatm*pi180)*r_earth |
---|
88 | if (cosfactp.lt.cosfactm) then |
---|
89 | hzone=sqrt(r_earth**2-cosfactp**2)-sqrt(r_earth**2-cosfactm**2) |
---|
90 | else |
---|
91 | hzone=sqrt(r_earth**2-cosfactm**2)-sqrt(r_earth**2-cosfactp**2) |
---|
92 | endif |
---|
93 | endif |
---|
94 | ! |
---|
95 | !C Surface are of a grid cell at a latitude ylat |
---|
96 | !*********************************************** |
---|
97 | ! |
---|
98 | gridarea=2.*pi*r_earth*hzone*dxoutl/360. |
---|
99 | |
---|
100 | ! for FLEXPART_WRF, dx & dy are in meters, and no cos(lat) is needed |
---|
101 | ! ??? maybe should incorporate map factor here, |
---|
102 | ! and also for areaeast & areanorth ??? |
---|
103 | ! gridarea=dxout*dyout |
---|
104 | |
---|
105 | do ix=0,numxgrid-1 |
---|
106 | area(ix,jy)=gridarea |
---|
107 | ! Volume = area x box height |
---|
108 | !*************************** |
---|
109 | |
---|
110 | volume(ix,jy,1)=area(ix,jy)*outheight(1) |
---|
111 | |
---|
112 | ! for FLEXPART_WRF, dx & dy are in meters, and no cos(lat) is needed |
---|
113 | areaeast(ix,jy,1)=dyoutl*r_earth*pi180*outheight(1) |
---|
114 | areanorth(ix,jy,1)=cos(ylat*pi180)*dxoutl*r_earth*pi180*outheight(1) |
---|
115 | ! areaeast(ix,jy,1)=dyout*outheight(1) |
---|
116 | ! areanorth(ix,jy,1)=dxout*outheight(1) |
---|
117 | |
---|
118 | do kz=2,numzgrid |
---|
119 | |
---|
120 | areaeast(ix,jy,kz)=dyoutl*r_earth*pi180*(outheight(kz)-outheight(kz-1)) |
---|
121 | areanorth(ix,jy,kz)=cos(ylat*pi180)*dxoutl*r_earth*pi180* & |
---|
122 | (outheight(kz)-outheight(kz-1)) |
---|
123 | ! areaeast(ix,jy,kz)=dyout*(outheight(kz)-outheight(kz-1)) |
---|
124 | ! areanorth(ix,jy,kz)=dxout*(outheight(kz)-outheight(kz-1)) |
---|
125 | |
---|
126 | volume(ix,jy,kz)=area(ix,jy)*(outheight(kz)-outheight(kz-1)) |
---|
127 | end do |
---|
128 | end do |
---|
129 | end do |
---|
130 | |
---|
131 | |
---|
132 | !****************************************************************** |
---|
133 | ! Determine average height of model topography in output grid cells |
---|
134 | !****************************************************************** |
---|
135 | |
---|
136 | ! Loop over all output grid cells |
---|
137 | !******************************** |
---|
138 | |
---|
139 | do jjy=0,numygrid-1 |
---|
140 | do iix=0,numxgrid-1 |
---|
141 | oroh=0. |
---|
142 | |
---|
143 | ! Take 100 samples of the topography in every grid cell |
---|
144 | !****************************************************** |
---|
145 | |
---|
146 | do j1=1,10 |
---|
147 | ! for FLEXPART_WRF, x & y coords are in meters, |
---|
148 | ! and the lon & lat variables below are in meters. |
---|
149 | ylat=outlat0+(float(jjy)+float(j1)/10.-0.05)*dyoutl !in degrees |
---|
150 | ! ymet=out_ym0+(real(jjy)+real(j1)/10.-0.05)*dyout |
---|
151 | ! yl=(ymet-ymet0)/dy |
---|
152 | do i1=1,10 |
---|
153 | xlon=outlon0+(float(iix)+float(i1)/10.-0.05)*dxoutl !in degrees |
---|
154 | call ll_to_xymeter_wrf(xlon,ylat,xmet,ymet) |
---|
155 | yl=(ymet-ymet0)/dy |
---|
156 | xl=(xmet-xmet0)/dx |
---|
157 | print*,'orotest',jjy,iix,j1,i1,xlon,xmet,xl,dx,dxoutl,outlon0 |
---|
158 | |
---|
159 | ! xmet=out_xm0+(real(iix)+real(i1)/10.-0.05)*dxout |
---|
160 | ! xl=(xmet-xmet0)/dx |
---|
161 | |
---|
162 | ! Determine the nest we are in |
---|
163 | !***************************** |
---|
164 | |
---|
165 | ngrid=0 |
---|
166 | do j=numbnests,1,-1 |
---|
167 | if ((xl.gt.xln(j)).and.(xl.lt.xrn(j)).and. & |
---|
168 | (yl.gt.yln(j)).and.(yl.lt.yrn(j))) then |
---|
169 | ngrid=j |
---|
170 | goto 43 |
---|
171 | endif |
---|
172 | end do |
---|
173 | 43 continue |
---|
174 | |
---|
175 | ! Determine (nested) grid coordinates and auxiliary parameters |
---|
176 | ! used for interpolation |
---|
177 | !************************************************************** |
---|
178 | |
---|
179 | if (ngrid.gt.0) then |
---|
180 | xtn=(xl-xln(ngrid))*xresoln(ngrid) |
---|
181 | ytn=(yl-yln(ngrid))*yresoln(ngrid) |
---|
182 | ix=int(xtn) |
---|
183 | jy=int(ytn) |
---|
184 | ddy=ytn-real(jy) |
---|
185 | ddx=xtn-real(ix) |
---|
186 | else |
---|
187 | ix=int(xl) |
---|
188 | jy=int(yl) |
---|
189 | ddy=yl-real(jy) |
---|
190 | ddx=xl-real(ix) |
---|
191 | endif |
---|
192 | ixp=ix+1 |
---|
193 | jyp=jy+1 |
---|
194 | rddx=1.-ddx |
---|
195 | rddy=1.-ddy |
---|
196 | p1=rddx*rddy |
---|
197 | p2=ddx*rddy |
---|
198 | p3=rddx*ddy |
---|
199 | p4=ddx*ddy |
---|
200 | |
---|
201 | if (ngrid.gt.0) then |
---|
202 | oroh=oroh+ p1*oron(ix ,jy ,ngrid) & |
---|
203 | + p2*oron(ixp,jy ,ngrid) & |
---|
204 | + p3*oron(ix ,jyp,ngrid) & |
---|
205 | + p4*oron(ixp,jyp,ngrid) |
---|
206 | else |
---|
207 | oroh=oroh+ p1*oro(ix ,jy) & |
---|
208 | + p2*oro(ixp,jy) & |
---|
209 | + p3*oro(ix ,jyp) & |
---|
210 | + p4*oro(ixp,jyp) |
---|
211 | endif |
---|
212 | end do |
---|
213 | end do |
---|
214 | |
---|
215 | ! Divide by the number of samples taken |
---|
216 | !************************************** |
---|
217 | |
---|
218 | oroout(iix,jjy)=oroh/100. |
---|
219 | |
---|
220 | end do |
---|
221 | end do |
---|
222 | |
---|
223 | |
---|
224 | ! if necessary allocate flux fields |
---|
225 | if (iflux.eq.1) then |
---|
226 | allocate(flux(6,0:numxgrid-1,0:numygrid-1,numzgrid, & |
---|
227 | 1:nspec,1:maxpointspec_act,1:nageclass),stat=stat) |
---|
228 | if (stat.ne.0) write(*,*)'ERROR: could not allocate flux array ' |
---|
229 | endif |
---|
230 | |
---|
231 | !write (*,*) 'allocating: in a sec',OHREA |
---|
232 | if (OHREA) then |
---|
233 | ! write (*,*) 'allocating: ',maxxOH,maxyOH,maxzOH |
---|
234 | allocate(OH_field(12,0:maxxOH-1,0:maxyOH-1,maxzOH),stat=stat) |
---|
235 | if (stat.ne.0) write(*,*)'ERROR: could not allocate OH array ' |
---|
236 | allocate(OH_field_height(7),stat=stat) |
---|
237 | if (stat.ne.0) write(*,*)'ERROR: could not allocate OH array ' |
---|
238 | endif |
---|
239 | ! gridunc,griduncn uncertainty of outputted concentrations |
---|
240 | allocate(gridunc(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, & |
---|
241 | maxpointspec_act,nclassunc,maxageclass),stat=stat) |
---|
242 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
243 | |
---|
244 | print*,'alloc gridunc',numxgrid-1,numygrid-1,numzgrid,maxspec, & |
---|
245 | maxpointspec_act,nclassunc,maxageclass |
---|
246 | if (ldirect.gt.0) then |
---|
247 | allocate(wetgridunc(0:numxgrid-1,0:numygrid-1,maxspec, & |
---|
248 | maxpointspec_act,nclassunc,maxageclass),stat=stat) |
---|
249 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
250 | allocate(drygridunc(0:numxgrid-1,0:numygrid-1,maxspec, & |
---|
251 | maxpointspec_act,nclassunc,maxageclass),stat=stat) |
---|
252 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
253 | allocate(drygridunc2(0:numxgrid-1,0:numygrid-1,maxspec, & |
---|
254 | maxpointspec_act,nclassunc,maxageclass),stat=stat) |
---|
255 | endif |
---|
256 | !write (*,*) 'Dimensions for fields', numxgrid,numygrid, & |
---|
257 | ! maxspec,maxpointspec_act,nclassunc,maxageclass |
---|
258 | |
---|
259 | write (*,*) ' Allocating fields for nested and global output (x,y): ', & |
---|
260 | max(numxgrid,numxgridn),max(numygrid,numygridn) |
---|
261 | |
---|
262 | ! allocate fields for concoutput with maximum dimension of outgrid |
---|
263 | ! and outgrid_nest |
---|
264 | allocate(gridsigma(0:max(numxgrid,numxgridn)-1, & |
---|
265 | 0:max(numygrid,numygridn)-1,numzgrid),stat=stat) |
---|
266 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
267 | allocate(grid(0:max(numxgrid,numxgridn)-1, & |
---|
268 | 0:max(numygrid,numygridn)-1,numzgrid),stat=stat) |
---|
269 | allocate(grid2(0:max(numxgrid,numxgridn)-1, 0:max(numygrid,numygridn)-1, & |
---|
270 | numzgrid,maxpointspec_act,nageclass),stat=stat) |
---|
271 | allocate(grid3(0:max(numxgrid,numxgridn)-1, 0:max(numygrid,numygridn)-1, & |
---|
272 | numzgrid,maxpointspec_act,nageclass),stat=stat) |
---|
273 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
274 | allocate(densityoutgrid(0:max(numxgrid,numxgridn)-1, & |
---|
275 | 0:max(numygrid,numygridn)-1,numzgrid),stat=stat) |
---|
276 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
277 | |
---|
278 | allocate(factor3d(0:max(numxgrid,numxgridn)-1, & |
---|
279 | 0:max(numygrid,numygridn)-1,numzgrid),stat=stat) |
---|
280 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
281 | allocate(sparse_dump_r(max(numxgrid,numxgridn)* & |
---|
282 | max(numygrid,numygridn)*numzgrid),stat=stat) |
---|
283 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
284 | allocate(sparse_dump_i(max(numxgrid,numxgridn)* & |
---|
285 | max(numygrid,numygridn)*numzgrid),stat=stat) |
---|
286 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
287 | |
---|
288 | ! deposition fields are only allocated for forward runs |
---|
289 | if (ldirect.gt.0) then |
---|
290 | allocate(wetgridsigma(0:max(numxgrid,numxgridn)-1, & |
---|
291 | 0:max(numygrid,numygridn)-1),stat=stat) |
---|
292 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
293 | allocate(drygridsigma(0:max(numxgrid,numxgridn)-1, & |
---|
294 | 0:max(numygrid,numygridn)-1),stat=stat) |
---|
295 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
296 | allocate(wetgrid(0:max(numxgrid,numxgridn)-1, & |
---|
297 | 0:max(numygrid,numygridn)-1),stat=stat) |
---|
298 | allocate(wetgrid2(0:max(numxgrid,numxgridn)-1, & |
---|
299 | 0:max(numygrid,numygridn)-1,maxpointspec_act,nageclass),stat=stat) |
---|
300 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
301 | allocate(drygrid(0:max(numxgrid,numxgridn)-1, & |
---|
302 | 0:max(numygrid,numygridn)-1),stat=stat) |
---|
303 | allocate(drygrid2(0:max(numxgrid,numxgridn)-1, & |
---|
304 | 0:max(numygrid,numygridn)-1,maxpointspec_act,nageclass),stat=stat) |
---|
305 | if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc' |
---|
306 | |
---|
307 | endif |
---|
308 | |
---|
309 | ! Initial condition field |
---|
310 | |
---|
311 | if (linit_cond.gt.0) then |
---|
312 | allocate(init_cond(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, & |
---|
313 | maxpointspec_act),stat=stat) |
---|
314 | if (stat.ne.0) write(*,*)'ERROR: could not allocate init_cond' |
---|
315 | endif |
---|
316 | |
---|
317 | !************************ |
---|
318 | ! Initialize output grids |
---|
319 | !************************ |
---|
320 | |
---|
321 | do ks=1,nspec |
---|
322 | do kp=1,maxpointspec_act |
---|
323 | do i=1,numreceptor |
---|
324 | ! Receptor points |
---|
325 | creceptor(i,ks)=0. |
---|
326 | end do |
---|
327 | do nage=1,nageclass |
---|
328 | do jy=0,numygrid-1 |
---|
329 | do ix=0,numxgrid-1 |
---|
330 | do l=1,nclassunc |
---|
331 | ! Deposition fields |
---|
332 | if (ldirect.gt.0) then |
---|
333 | wetgridunc(ix,jy,ks,kp,l,nage)=0. |
---|
334 | drygridunc(ix,jy,ks,kp,l,nage)=0. |
---|
335 | endif |
---|
336 | do kz=1,numzgrid |
---|
337 | if (iflux.eq.1) then |
---|
338 | ! Flux fields |
---|
339 | do i=1,5 |
---|
340 | flux(i,ix,jy,kz,ks,kp,nage)=0. |
---|
341 | end do |
---|
342 | endif |
---|
343 | ! Initial condition field |
---|
344 | if ((l.eq.1).and.(nage.eq.1).and.(linit_cond.gt.0)) & |
---|
345 | init_cond(ix,jy,kz,ks,kp)=0. |
---|
346 | ! Concentration fields |
---|
347 | gridunc(ix,jy,kz,ks,kp,l,nage)=0. |
---|
348 | end do |
---|
349 | end do |
---|
350 | end do |
---|
351 | end do |
---|
352 | end do |
---|
353 | end do |
---|
354 | end do |
---|
355 | |
---|
356 | end subroutine outgrid_init_reg |
---|