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 | |
---|
22 | program flexpart |
---|
23 | |
---|
24 | !***************************************************************************** |
---|
25 | ! * |
---|
26 | ! This is the Lagrangian Particle Dispersion Model FLEXPART. * |
---|
27 | ! The main program manages the reading of model run specifications, etc. * |
---|
28 | ! All actual computing is done within subroutine timemanager. * |
---|
29 | ! * |
---|
30 | ! Author: A. Stohl * |
---|
31 | ! * |
---|
32 | ! 18 May 1996 * |
---|
33 | ! * |
---|
34 | !***************************************************************************** |
---|
35 | ! * |
---|
36 | ! Variables: * |
---|
37 | ! * |
---|
38 | ! Constants: * |
---|
39 | ! * |
---|
40 | !***************************************************************************** |
---|
41 | |
---|
42 | use point_mod |
---|
43 | use par_mod |
---|
44 | use com_mod |
---|
45 | use conv_mod |
---|
46 | use random_mod, only: gasdev1 |
---|
47 | #if (defined _OPENMP) |
---|
48 | use omp_lib, only: OMP_GET_MAX_THREADS |
---|
49 | #endif |
---|
50 | |
---|
51 | #ifdef NETCDF_OUTPUT |
---|
52 | use nc_output_mod, only: writeheader_ncdf |
---|
53 | #endif |
---|
54 | |
---|
55 | |
---|
56 | implicit none |
---|
57 | |
---|
58 | integer :: i,j,ix,jy,inest |
---|
59 | integer :: idummy = -320 |
---|
60 | #if (defined _OPENMP) |
---|
61 | integer :: numthreads |
---|
62 | #endif |
---|
63 | character(256) :: pathfile |
---|
64 | integer :: clck_counts_beg, clck_counts_end, clck_rate, tins |
---|
65 | |
---|
66 | ! for measuring wall time |
---|
67 | call system_clock (clck_counts_beg,clck_rate) |
---|
68 | |
---|
69 | ! Generate a large number of random numbers |
---|
70 | !****************************************** |
---|
71 | |
---|
72 | do i=1,maxrand-1,2 |
---|
73 | call gasdev1(idummy,rannumb(i),rannumb(i+1)) |
---|
74 | end do |
---|
75 | call gasdev1(idummy,rannumb(maxrand),rannumb(maxrand-1)) |
---|
76 | |
---|
77 | ! Print the GPL License statement |
---|
78 | !******************************************************* |
---|
79 | print*,'Welcome to FLEXPART Version 9.1 (Build 20130723)' |
---|
80 | print*,'FLEXPART is free software released under the GNU General Public License.' |
---|
81 | |
---|
82 | ! Read the pathnames where input/output files are stored |
---|
83 | !******************************************************* |
---|
84 | |
---|
85 | select case (iargc()) |
---|
86 | case (1) |
---|
87 | call getarg(1,pathfile) |
---|
88 | case (0) |
---|
89 | write(pathfile,'(a11)') './pathnames' |
---|
90 | end select |
---|
91 | call readpaths(pathfile) |
---|
92 | |
---|
93 | ! Read the user specifications for the current model run |
---|
94 | !******************************************************* |
---|
95 | |
---|
96 | call readcommand |
---|
97 | |
---|
98 | ! Check how many threads are available in parallel section |
---|
99 | !******************************************************** |
---|
100 | |
---|
101 | #if (defined _OPENMP) |
---|
102 | numthreads = OMP_GET_MAX_THREADS() |
---|
103 | if (numthreads.gt.1) then |
---|
104 | write(*,*) |
---|
105 | write(*,*) "******************** INFO **************************" |
---|
106 | write(*,*) "* FLEXPART running in parallel mode *" |
---|
107 | write(*,*) "* Number of threads:", numthreads, ". *" |
---|
108 | write(*,*) "******************************************************" |
---|
109 | write(*,*) |
---|
110 | endif |
---|
111 | #endif |
---|
112 | |
---|
113 | |
---|
114 | ! Read the age classes to be used |
---|
115 | !******************************** |
---|
116 | |
---|
117 | call readageclasses |
---|
118 | |
---|
119 | |
---|
120 | ! Read, which wind fields are available within the modelling period |
---|
121 | !****************************************************************** |
---|
122 | |
---|
123 | call readavailable |
---|
124 | |
---|
125 | |
---|
126 | ! Read the model grid specifications, |
---|
127 | ! both for the mother domain and eventual nests |
---|
128 | !********************************************** |
---|
129 | |
---|
130 | call gridcheck |
---|
131 | call gridcheck_nests |
---|
132 | |
---|
133 | |
---|
134 | ! Read the output grid specifications |
---|
135 | !************************************ |
---|
136 | |
---|
137 | call readoutgrid |
---|
138 | if (nested_output.eq.1) call readoutgrid_nest |
---|
139 | |
---|
140 | |
---|
141 | ! Read the receptor points for which extra concentrations are to be calculated |
---|
142 | !***************************************************************************** |
---|
143 | |
---|
144 | call readreceptors |
---|
145 | |
---|
146 | |
---|
147 | ! Read the physico-chemical species property table |
---|
148 | !************************************************* |
---|
149 | !SEC: now only needed SPECIES are read in readreleases.f |
---|
150 | !call readspecies |
---|
151 | |
---|
152 | |
---|
153 | ! Read the landuse inventory |
---|
154 | !*************************** |
---|
155 | |
---|
156 | call readlanduse |
---|
157 | |
---|
158 | |
---|
159 | ! Assign fractional cover of landuse classes to each ECMWF grid point |
---|
160 | !******************************************************************** |
---|
161 | |
---|
162 | call assignland |
---|
163 | |
---|
164 | |
---|
165 | |
---|
166 | ! Read the coordinates of the release locations |
---|
167 | !********************************************** |
---|
168 | |
---|
169 | call readreleases |
---|
170 | |
---|
171 | ! Read and compute surface resistances to dry deposition of gases |
---|
172 | !**************************************************************** |
---|
173 | |
---|
174 | call readdepo |
---|
175 | |
---|
176 | |
---|
177 | ! Convert the release point coordinates from geografical to grid coordinates |
---|
178 | !*************************************************************************** |
---|
179 | |
---|
180 | call coordtrafo |
---|
181 | |
---|
182 | |
---|
183 | ! Initialize all particles to non-existent |
---|
184 | !***************************************** |
---|
185 | |
---|
186 | do j=1,maxpart |
---|
187 | itra1(j)=-999999999 |
---|
188 | end do |
---|
189 | |
---|
190 | ! For continuation of previous run, read in particle positions |
---|
191 | !************************************************************* |
---|
192 | |
---|
193 | if (ipin.eq.1) then |
---|
194 | call readpartpositions |
---|
195 | else |
---|
196 | numpart=0 |
---|
197 | numparticlecount=0 |
---|
198 | endif |
---|
199 | |
---|
200 | |
---|
201 | ! Calculate volume, surface area, etc., of all output grid cells |
---|
202 | ! Allocate fluxes and OHfield if necessary |
---|
203 | !*************************************************************** |
---|
204 | |
---|
205 | call outgrid_init |
---|
206 | if (nested_output.eq.1) call outgrid_init_nest |
---|
207 | |
---|
208 | |
---|
209 | ! Read the OH field |
---|
210 | !****************** |
---|
211 | |
---|
212 | if (OHREA.eqv..TRUE.) & |
---|
213 | call readOHfield |
---|
214 | |
---|
215 | ! Write basic information on the simulation to a file "header" |
---|
216 | ! and open files that are to be kept open throughout the simulation |
---|
217 | !****************************************************************** |
---|
218 | |
---|
219 | if (lncdfout) then |
---|
220 | #ifdef NETCDF_OUTPUT |
---|
221 | call writeheader_ncdf(lnest = .false.) |
---|
222 | #endif |
---|
223 | else |
---|
224 | call writeheader |
---|
225 | end if |
---|
226 | |
---|
227 | if (nested_output.eq.1) then |
---|
228 | if (lncdfout) then |
---|
229 | #ifdef NETCDF_OUTPUT |
---|
230 | call writeheader_ncdf(lnest = .true.) |
---|
231 | #endif |
---|
232 | else |
---|
233 | call writeheader_nest |
---|
234 | endif |
---|
235 | endif |
---|
236 | |
---|
237 | if (.NOT.lncdfout) then |
---|
238 | open(unitdates,file=path(2)(1:length(2))//'dates') |
---|
239 | end if |
---|
240 | |
---|
241 | call openreceptors |
---|
242 | if ((iout.eq.4).or.(iout.eq.5)) call openouttraj |
---|
243 | |
---|
244 | |
---|
245 | ! Releases can only start and end at discrete times (multiples of lsynctime) |
---|
246 | !*************************************************************************** |
---|
247 | |
---|
248 | do i=1,numpoint |
---|
249 | ireleasestart(i)=nint(real(ireleasestart(i))/ & |
---|
250 | real(lsynctime))*lsynctime |
---|
251 | ireleaseend(i)=nint(real(ireleaseend(i))/ & |
---|
252 | real(lsynctime))*lsynctime |
---|
253 | end do |
---|
254 | |
---|
255 | |
---|
256 | ! Initialize cloud-base mass fluxes for the convection scheme |
---|
257 | !************************************************************ |
---|
258 | |
---|
259 | do jy=0,nymin1 |
---|
260 | do ix=0,nxmin1 |
---|
261 | cbaseflux(ix,jy)=0. |
---|
262 | end do |
---|
263 | end do |
---|
264 | do inest=1,numbnests |
---|
265 | do jy=0,nyn(inest)-1 |
---|
266 | do ix=0,nxn(inest)-1 |
---|
267 | cbasefluxn(ix,jy,inest)=0. |
---|
268 | end do |
---|
269 | end do |
---|
270 | end do |
---|
271 | |
---|
272 | |
---|
273 | ! Calculate particle trajectories |
---|
274 | !******************************** |
---|
275 | |
---|
276 | write(*,*) " Starting calculations" |
---|
277 | call timemanager |
---|
278 | |
---|
279 | |
---|
280 | write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLEXPART MODEL RUN!' |
---|
281 | write(*,*) |
---|
282 | |
---|
283 | ! output wall time |
---|
284 | call system_clock(clck_counts_end,clck_rate) |
---|
285 | tins=(clck_counts_end - clck_counts_beg)/real(clck_rate) |
---|
286 | print*,'Wall time ',tins,'s, ',tins/60,'min, ',tins/3600,'h.' |
---|
287 | |
---|
288 | end program flexpart |
---|