CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
FLEXPART.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 
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 
47  implicit none
48 
49  integer :: i,j,ix,jy,inest
50  integer :: idummy = -320
51 
52  integer :: metdata_format = unknown_metdata
53 
54 #ifdef TESTUSEGETFPFIELDS
55  ! Data structures used to pirate the data structures filled from
56  ! pathnames and AVAILABLE and point to test .fp binary files.
57 
58  integer :: idx_wf ! Stores the index number of the met/.fp file
59 
60  ! Full path to the .fp files location
61  !character(len=*), parameter :: FPDIR = '/home/morton/basic_ec_nc_testing/'
62  character(len=*), parameter :: fpdir = '/home/morton/'
63  !character(len=*), parameter :: FPDIR = '/home/morton/FPForward/'
64  !character(len=*), parameter :: FPDIR = '/home/morton/FPTinyFwd/'
65  !character(len=*), parameter :: FPDIR = '/home/morton/ecmwf_tiny_fp/'
66 
67  character(len=512) :: fpname ! Stores name of a single .fp file
68 #endif
69 
70  ! Generate a large number of random numbers
71  !******************************************
72 
73  do i=1,maxrand-1,2
74  call gasdev1(idummy,rannumb(i),rannumb(i+1))
75  end do
76  call gasdev1(idummy,rannumb(maxrand),rannumb(maxrand-1))
77 
78  ! Print the GPL License statement
79  !*******************************************************
80 #if defined WITH_CTBTO_PATCHES
81  print*,'Welcome to FLEXPART Version 9.0 CTBTO'
82 #else
83  print*,'Welcome to FLEXPART Version 9.0'
84 #endif
85 
86  print*,'FLEXPART is free software released under the GNU Genera'// &
87  'l Public License.'
88 
89  ! Read the pathnames where input/output files are stored
90  !*******************************************************
91 
92  call readpaths
93 
94  ! Read the user specifications for the current model run
95  !*******************************************************
96 
97  call readcommand
98 
99 
100  ! Read the age classes to be used
101  !********************************
102 
103  call readageclasses
104 
105 
106  ! Read, which wind fields are available within the modelling period
107  !******************************************************************
108 
109  call readavailable
110 
111 
112 
113 
114 
115 #ifdef TESTUSEGETFPFIELDS
116  ! THIS IS ONLY USED FOR TESTING
117  ! Hard code the path to .fp binary files into the data structures as if they
118  ! were read in from pathnames and the AVAILABLE file up above in
119  ! readpaths() and readavailable(). In short, this replaces the path and
120  ! names of the GRIB met files with those of the .fp binary files.
121  ! THIS IS ONLY USED FOR TESTING
122 
123  ! I think it might be problematic that the AVAILABLE file just might have
124  ! more files than there are .fp files, because the generation of the .fp
125  ! files may have taken place with a subset of the full AVAILABLE file list.
126 
127  ! These variables, used below, come from com_mod
128  !
129  ! --- numbwf, wfname, path, length
130 
131  ! This test code needs to come after gridcheck() and gridcheck_nests(), since
132  ! at the time those are called, the data structures still need to hold path
133  ! names to GRIB files
134 
135  print *, 'Before overwriting...'
136  print *, 'path(3): ', path(3)
137  print *, 'length(3): ', length(3)
138 
139  ! Get and store the base name of the metfiles read in from AVAILABLE
140  DO idx_wf = 1, numbwf
141  wfname(idx_wf) = trim(wfname(idx_wf))
142  print *, wfname(idx_wf)
143  END DO
144 
145  ! The values in path and length arrays were meant to store met file path
146  ! and length of path. Replace with the .fp file path and length of path
147  path(3) = fpdir
148  length(3) = len(trim(path(3)))
149 
150  print *, 'After overwriting...'
151  print *, 'path(3): ', path(3)
152  print *, 'length(3): ', length(3)
153 
154  ! Add the .fp extension to each of the met file base names
155  DO idx_wf = 1, numbwf
156  fpname = trim(wfname(idx_wf)) // ".fp"
157  wfname(idx_wf) = trim(fpname)
158  print *, wfname(idx_wf)
159  END DO
160 
161  ! At this point, the path and filename data structures have been overwritten
162  ! to point to the .fp files, rather than the metfiles
163 
164  !!!!!!!!!!!!!!! CALL TO fpgridcheck()
165  CALL fpgridcheck
166 #else
167 
168  if ( preprocessed_metdata.eq.1 ) then
169  call fpgridcheck
170  else
171  ! Detect metdata format
172  !******************************************************************
173  call detectformat(metdata_format)
174 
175  if (metdata_format.eq.ecmwf_metdata) then
176  print*,'ECMWF metdata detected'
177  elseif (metdata_format.eq.gfs_metdata) then
178  print*,'NCEP metdata detected'
179  else
180  stop 'Unknown metdata format'
181  endif
182 
183  ! Read the model grid specifications,
184  ! both for the mother domain and eventual nests
185  !**********************************************
186  if (metdata_format.eq.ecmwf_metdata) call gridcheck_ecmwf
187  if (metdata_format.eq.gfs_metdata) call gridcheck_gfs
188  call gridcheck_nests
189  endif
190 #endif
191 
192  ! Read the output grid specifications
193  !************************************
194  call readoutgrid
195  if (nested_output.eq.1) call readoutgrid_nest
196 
197 
198  ! Read the receptor points for which extra concentrations are to be calculated
199  !*****************************************************************************
200  call readreceptors
201 
202 
203  ! Read the physico-chemical species property table
204  !*************************************************
205  !SEC: now only needed SPECIES are read in readreleases.f
206  !call readspecies
207 
208 
209  ! Read the landuse inventory
210  !***************************
211  call readlanduse
212 
213 
214  ! Assign fractional cover of landuse classes to each ECMWF grid point
215  !********************************************************************
216  call assignland
217 
218 
219 
220  ! Read the coordinates of the release locations
221  !**********************************************
222 
223  call readreleases
224 
225  ! Read and compute surface resistances to dry deposition of gases
226  !****************************************************************
227 
228  call readdepo
229 
230 
231  ! Convert the release point coordinates from geografical to grid coordinates
232  !***************************************************************************
233 
234  call coordtrafo
235 
236 
237  ! Initialize all particles to non-existent
238  !*****************************************
239 
240  do j=1,maxpart
241  itra1(j)=-999999999
242  end do
243 
244  ! For continuation of previous run, read in particle positions
245  !*************************************************************
246 
247  if (ipin.eq.1) then
248  call readpartpositions
249  else
250  numpart=0
251  numparticlecount=0
252  endif
253 
254 
255  ! Calculate volume, surface area, etc., of all output grid cells
256  ! Allocate fluxes and OHfield if necessary
257  !***************************************************************
258 
259  call outgrid_init
260  if (nested_output.eq.1) call outgrid_init_nest
261 
262 
263  ! Read the OH field
264  !******************
265 
266  if (ohrea.eqv..true.) &
267  call readohfield
268 
269  ! Write basic information on the simulation to a file "header"
270  ! and open files that are to be kept open throughout the simulation
271  !******************************************************************
272 
273  call writeheader
274  if (nested_output.eq.1) call writeheader_nest
275  open(unitdates,file=path(2)(1:length(2))//'dates')
276  call openreceptors
277  if ((iout.eq.4).or.(iout.eq.5)) call openouttraj
278 
279 
280  ! Releases can only start and end at discrete times (multiples of lsynctime)
281  !***************************************************************************
282 
283  do i=1,numpoint
284  ireleasestart(i)=nint(real(ireleasestart(i))/ &
285  real(lsynctime))*lsynctime
286  ireleaseend(i)=nint(real(ireleaseend(i))/ &
287  real(lsynctime))*lsynctime
288  end do
289 
290 
291  ! Initialize cloud-base mass fluxes for the convection scheme
292  !************************************************************
293 
294  do jy=0,nymin1
295  do ix=0,nxmin1
296  cbaseflux(ix,jy)=0.
297  end do
298  end do
299  do inest=1,numbnests
300  do jy=0,nyn(inest)-1
301  do ix=0,nxn(inest)-1
302  cbasefluxn(ix,jy,inest)=0.
303  end do
304  end do
305  end do
306 
307 
308  ! Calculate particle trajectories
309  !********************************
310 
311 
312  call timemanager(metdata_format)
313 
314 
315  write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
316  &XPART MODEL RUN!'
317 
318 end program flexpart
subroutine readpaths
Definition: readpaths.f90:22
subroutine assignland
Definition: assignland.f90:22
subroutine readdepo
Definition: readdepo.f90:22
subroutine readreceptors
subroutine openouttraj
Definition: openouttraj.f90:22
subroutine readoutgrid
Definition: readoutgrid.f90:22
subroutine writeheader_nest
subroutine writeheader
Definition: writeheader.f90:22
subroutine readavailable
subroutine coordtrafo
Definition: coordtrafo.f90:22
subroutine readcommand
Definition: readcommand.F90:22
subroutine readageclasses
subroutine gridcheck_ecmwf
Definition: gridcheck.F90:22
subroutine readreleases
subroutine timemanager(metdata_format)
Definition: timemanager.F90:22
subroutine gridcheck_nests
subroutine gridcheck_gfs
subroutine gasdev1(idum, random1, random2)
Definition: random.f90:83
subroutine readoutgrid_nest
subroutine openreceptors
subroutine detectformat(metdata_format)
program flexpart
Definition: FLEXPART.F90:22
subroutine outgrid_init_nest
subroutine readpartpositions
subroutine readlanduse
Definition: readlanduse.f90:22
subroutine readohfield
Definition: readOHfield.f90:22
subroutine outgrid_init
subroutine fpgridcheck
Definition: fpgridcheck.F90:2