source: branches/FLEXPART_9.1.3/src/FLEXPART.f90 @ 13

Last change on this file since 13 was 13, checked in by saeck, 11 years ago

update to wetdepo.f90

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