source: flexpart.git/src/FLEXPART.f90 @ e200b7a

10.4.1_peseiFPv9.3.1FPv9.3.1b_testingFPv9.3.2GFS_025NetCDFbugfixes+enhancementsdepositiondevflexpart-noresmfp9.3.1-20161214-nc4grib2nc4_repairinputlistlaptoprelease-10release-10.4.1scaling-bugsvn-petrasvn-trunkunivie
Last change on this file since e200b7a was e200b7a, checked in by Matthias Langer <matthias.langer@…>, 11 years ago

git-svn-id: http://flexpart.flexpart.eu:8088/svn/FlexPart90/trunk@4 ef8cc7e1-21b7-489e-abab-c1baa636049d

  • Property mode set to 100644
File size: 7.3 KB
RevLine 
[e200b7a]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
22program 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  ! Generate a large number of random numbers
53  !******************************************
54
55  do i=1,maxrand-1,2
56    call gasdev1(idummy,rannumb(i),rannumb(i+1))
57  end do
58  call gasdev1(idummy,rannumb(maxrand),rannumb(maxrand-1))
59
60  ! Print the GPL License statement
61  !*******************************************************
62  print*,'Welcome to FLEXPART Version 9.0'
63  print*,'FLEXPART is free software released under the GNU Genera'// &
64       'l Public License.'
65
66  ! Read the pathnames where input/output files are stored
67  !*******************************************************
68
69  call readpaths
70
71  ! Read the user specifications for the current model run
72  !*******************************************************
73
74  call readcommand
75
76
77  ! Read the age classes to be used
78  !********************************
79
80  call readageclasses
81
82
83  ! Read, which wind fields are available within the modelling period
84  !******************************************************************
85
86  call readavailable
87
88
89  ! Read the model grid specifications,
90  ! both for the mother domain and eventual nests
91  !**********************************************
92
93  call gridcheck
94  call gridcheck_nests
95
96
97  ! Read the output grid specifications
98  !************************************
99
100  call readoutgrid
101  if (nested_output.eq.1) call readoutgrid_nest
102
103
104  ! Read the receptor points for which extra concentrations are to be calculated
105  !*****************************************************************************
106
107  call readreceptors
108
109
110  ! Read the physico-chemical species property table
111  !*************************************************
112  !SEC: now only needed SPECIES are read in readreleases.f
113  !call readspecies
114
115
116  ! Read the landuse inventory
117  !***************************
118
119  call readlanduse
120
121
122  ! Assign fractional cover of landuse classes to each ECMWF grid point
123  !********************************************************************
124
125  call assignland
126
127
128
129  ! Read the coordinates of the release locations
130  !**********************************************
131
132  call readreleases
133
134  ! Read and compute surface resistances to dry deposition of gases
135  !****************************************************************
136
137  call readdepo
138
139
140  ! Convert the release point coordinates from geografical to grid coordinates
141  !***************************************************************************
142
143  call coordtrafo
144
145
146  ! Initialize all particles to non-existent
147  !*****************************************
148
149  do j=1,maxpart
150    itra1(j)=-999999999
151  end do
152
153  ! For continuation of previous run, read in particle positions
154  !*************************************************************
155
156  if (ipin.eq.1) then
157    call readpartpositions
158  else
159    numpart=0
160    numparticlecount=0
161  endif
162
163
164  ! Calculate volume, surface area, etc., of all output grid cells
165  ! Allocate fluxes and OHfield if necessary
166  !***************************************************************
167
168  call outgrid_init
169  if (nested_output.eq.1) call outgrid_init_nest
170
171
172  ! Read the OH field
173  !******************
174
175  if (OHREA.eqv..TRUE.) &
176       call readOHfield
177
178  ! Write basic information on the simulation to a file "header"
179  ! and open files that are to be kept open throughout the simulation
180  !******************************************************************
181
182  call writeheader
183  if (nested_output.eq.1) call writeheader_nest
184  open(unitdates,file=path(2)(1:length(2))//'dates')
185  call openreceptors
186  if ((iout.eq.4).or.(iout.eq.5)) call openouttraj
187
188
189  ! Releases can only start and end at discrete times (multiples of lsynctime)
190  !***************************************************************************
191
192  do i=1,numpoint
193    ireleasestart(i)=nint(real(ireleasestart(i))/ &
194         real(lsynctime))*lsynctime
195    ireleaseend(i)=nint(real(ireleaseend(i))/ &
196         real(lsynctime))*lsynctime
197  end do
198
199
200  ! Initialize cloud-base mass fluxes for the convection scheme
201  !************************************************************
202
203  do jy=0,nymin1
204    do ix=0,nxmin1
205      cbaseflux(ix,jy)=0.
206    end do
207  end do
208  do inest=1,numbnests
209    do jy=0,nyn(inest)-1
210      do ix=0,nxn(inest)-1
211        cbasefluxn(ix,jy,inest)=0.
212      end do
213    end do
214  end do
215
216
217  ! Calculate particle trajectories
218  !********************************
219
220  call timemanager
221
222
223  write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
224       &XPART MODEL RUN!'
225
226end program flexpart
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG