CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
outgrid_init_nest.f90
Go to the documentation of this file.
1 !**********************************************************************
2 ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 *
3 ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, *
4 ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann *
5 ! *
6 ! This file is part of FLEXPART. *
7 ! *
8 ! FLEXPART is free software: you can redistribute it and/or modify *
9 ! it under the terms of the GNU General Public License as published by*
10 ! the Free Software Foundation, either version 3 of the License, or *
11 ! (at your option) any later version. *
12 ! *
13 ! FLEXPART is distributed in the hope that it will be useful, *
14 ! but WITHOUT ANY WARRANTY; without even the implied warranty of *
15 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
16 ! GNU General Public License for more details. *
17 ! *
18 ! You should have received a copy of the GNU General Public License *
19 ! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
20 !**********************************************************************
21 
23 
24  !*****************************************************************************
25  ! *
26  ! This routine calculates, for each grid cell of the output nest, the *
27  ! volume and the surface area. *
28  ! *
29  ! Author: A. Stohl *
30  ! *
31  ! 30 August 2004 *
32  ! *
33  !*****************************************************************************
34  ! *
35  ! Variables: *
36  ! *
37  ! arean surface area of all output nest cells *
38  ! volumen volumes of all output nest cells *
39  ! *
40  !*****************************************************************************
41 
42  use unc_mod
43  use outg_mod
44  use par_mod
45  use com_mod
46 
47  implicit none
48 
49  integer :: ix,jy,kz,ks,kp,nage,l,iix,jjy,ixp,jyp,i1,j1,j,ngrid
50  integer :: stat
51  real :: ylat,gridarea,ylatp,ylatm,hzone,cosfactm,cosfactp
52  real :: xlon,xl,yl,ddx,ddy,rddx,rddy,p1,p2,p3,p4,xtn,ytn,oroh
53  real,parameter :: eps=nxmax/3.e5
54 
55 
56 
57  ! gridunc,griduncn uncertainty of outputted concentrations
58  allocate(griduncn(0:numxgridn-1,0:numygridn-1,numzgrid,maxspec, &
59  maxpointspec_act,nclassunc,maxageclass),stat=stat)
60  if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
61 
62  if (ldirect.gt.0) then
63  allocate(wetgriduncn(0:numxgridn-1,0:numygridn-1,maxspec, &
64  maxpointspec_act,nclassunc,maxageclass),stat=stat)
65  if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
66  allocate(drygriduncn(0:numxgridn-1,0:numygridn-1,maxspec, &
67  maxpointspec_act,nclassunc,maxageclass),stat=stat)
68  if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
69  endif
70 
71  ! Compute surface area and volume of each grid cell: area, volume;
72  ! and the areas of the northward and eastward facing walls: areaeast, areanorth
73  !***********************************************************************
74 
75  do jy=0,numygridn-1
76  ylat=outlat0n+(real(jy)+0.5)*dyoutn
77  ylatp=ylat+0.5*dyoutn
78  ylatm=ylat-0.5*dyoutn
79  if ((ylatm.lt.0).and.(ylatp.gt.0.)) then
80  hzone=dyoutn*r_earth*pi180
81  else
82 
83  ! Calculate area of grid cell with formula M=2*pi*R*h*dx/360,
84  ! see Netz, Formeln der Mathematik, 5. Auflage (1983), p.90
85  !************************************************************
86 
87  cosfactp=cos(ylatp*pi180)
88  cosfactm=cos(ylatm*pi180)
89  if (cosfactp.lt.cosfactm) then
90  hzone=sqrt(1-cosfactp**2)- &
91  sqrt(1-cosfactm**2)
92  hzone=hzone*r_earth
93  else
94  hzone=sqrt(1-cosfactm**2)- &
95  sqrt(1-cosfactp**2)
96  hzone=hzone*r_earth
97  endif
98  endif
99 
100 
101 
102  ! Surface are of a grid cell at a latitude ylat
103  !**********************************************
104 
105  gridarea=2.*pi*r_earth*hzone*dxoutn/360.
106 
107  do ix=0,numxgridn-1
108  arean(ix,jy)=gridarea
109 
110  ! Volume = area x box height
111  !***************************
112 
113  volumen(ix,jy,1)=arean(ix,jy)*outheight(1)
114  do kz=2,numzgrid
115  volumen(ix,jy,kz)=arean(ix,jy)*(outheight(kz)-outheight(kz-1))
116  end do
117  end do
118  end do
119 
120 
121  !**************************************************************************
122  ! Determine average height of model topography in nesteed output grid cells
123  !**************************************************************************
124 
125  ! Loop over all output grid cells
126  !********************************
127 
128  do jjy=0,numygridn-1
129  do iix=0,numxgridn-1
130  oroh=0.
131 
132  ! Take 100 samples of the topography in every grid cell
133  !******************************************************
134 
135  do j1=1,10
136  ylat=outlat0n+(real(jjy)+real(j1)/10.-0.05)*dyoutn
137  yl=(ylat-ylat0)/dy
138  do i1=1,10
139  xlon=outlon0n+(real(iix)+real(i1)/10.-0.05)*dxoutn
140  xl=(xlon-xlon0)/dx
141 
142  ! Determine the nest we are in
143  !*****************************
144 
145  ngrid=0
146  do j=numbnests,1,-1
147  if ((xl.gt.xln(j)+eps).and.(xl.lt.xrn(j)-eps).and. &
148  (yl.gt.yln(j)+eps).and.(yl.lt.yrn(j)-eps)) then
149  ngrid=j
150  goto 43
151  endif
152  end do
153 43 continue
154 
155  ! Determine (nested) grid coordinates and auxiliary parameters used for interpolation
156  !*****************************************************************************
157 
158  if (ngrid.gt.0) then
159  xtn=(xl-xln(ngrid))*xresoln(ngrid)
160  ytn=(yl-yln(ngrid))*yresoln(ngrid)
161  ix=int(xtn)
162  jy=int(ytn)
163  ddy=ytn-real(jy)
164  ddx=xtn-real(ix)
165  else
166  ix=int(xl)
167  jy=int(yl)
168  ddy=yl-real(jy)
169  ddx=xl-real(ix)
170  endif
171  ixp=ix+1
172  jyp=jy+1
173  rddx=1.-ddx
174  rddy=1.-ddy
175  p1=rddx*rddy
176  p2=ddx*rddy
177  p3=rddx*ddy
178  p4=ddx*ddy
179 
180  if (ngrid.gt.0) then
181  oroh=oroh+p1*oron(ix ,jy ,ngrid) &
182  + p2*oron(ixp,jy ,ngrid) &
183  + p3*oron(ix ,jyp,ngrid) &
184  + p4*oron(ixp,jyp,ngrid)
185  else
186  oroh=oroh+p1*oro(ix ,jy) &
187  + p2*oro(ixp,jy) &
188  + p3*oro(ix ,jyp) &
189  + p4*oro(ixp,jyp)
190  endif
191  end do
192  end do
193 
194  ! Divide by the number of samples taken
195  !**************************************
196 
197  orooutn(iix,jjy)=oroh/100.
198  end do
199  end do
200 
201 
202 
203  !*******************************
204  ! Initialization of output grids
205  !*******************************
206 
207  do kp=1,maxpointspec_act
208  do ks=1,nspec
209  do nage=1,nageclass
210  do jy=0,numygridn-1
211  do ix=0,numxgridn-1
212  do l=1,nclassunc
213  ! Deposition fields
214  if (ldirect.gt.0) then
215  wetgriduncn(ix,jy,ks,kp,l,nage)=0.
216  drygriduncn(ix,jy,ks,kp,l,nage)=0.
217  endif
218  ! Concentration fields
219  do kz=1,numzgrid
220  griduncn(ix,jy,kz,ks,kp,l,nage)=0.
221  end do
222  end do
223  end do
224  end do
225  end do
226  end do
227  end do
228 
229 
230 end subroutine outgrid_init_nest
subroutine outgrid_init_nest