source: flexpart.git/src/FLEXPART.f90 @ 7999df47

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 7999df47 was db712a8, checked in by Espen Sollum ATMOS <eso@…>, 8 years ago

Completed handling of nested wind fields with cloud water (for wet deposition).

  • Property mode set to 100644
File size: 13.0 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  use netcdf_output_mod, only: writeheader_netcdf
47  use random_mod, only: gasdev1
48
49  implicit none
50
51  integer :: i,j,ix,jy,inest
52  integer :: idummy = -320
53  character(len=256) :: inline_options  !pathfile, flexversion, arg2
54
55
56  ! Initialize arrays in com_mod
57  !*****************************
58  call com_mod_allocate_part(maxpart)
59 
60
61  ! Generate a large number of random numbers
62  !******************************************
63
64  do i=1,maxrand-1,2
65    call gasdev1(idummy,rannumb(i),rannumb(i+1))
66  end do
67  call gasdev1(idummy,rannumb(maxrand),rannumb(maxrand-1))
68
69  ! FLEXPART version string
70  flexversion_major = '10' ! Major version number, also used for species file names
71  flexversion='Version '//trim(flexversion_major)//'.0beta (2015-05-01)'
72  verbosity=0
73
74  ! Read the pathnames where input/output files are stored
75  !*******************************************************
76
77  inline_options='none'
78  select case (iargc())
79  case (2)
80    call getarg(1,arg1)
81    pathfile=arg1
82    call getarg(2,arg2)
83    inline_options=arg2
84  case (1)
85    call getarg(1,arg1)
86    pathfile=arg1
87    if (arg1(1:1).eq.'-') then
88      write(pathfile,'(a11)') './pathnames'
89      inline_options=arg1
90    endif
91  case (0)
92    write(pathfile,'(a11)') './pathnames'
93  end select
94 
95  ! Print the GPL License statement
96  !*******************************************************
97  print*,'Welcome to FLEXPART ', trim(flexversion)
98  print*,'FLEXPART is free software released under the GNU General Public License.'
99 
100  if (inline_options(1:1).eq.'-') then
101    if (trim(inline_options).eq.'-v'.or.trim(inline_options).eq.'-v1') then
102       print*, 'Verbose mode 1: display detailed information during run'
103       verbosity=1
104    endif
105    if (trim(inline_options).eq.'-v2') then
106       print*, 'Verbose mode 2: display more detailed information during run'
107       verbosity=2
108    endif
109    if (trim(inline_options).eq.'-i') then
110       print*, 'Info mode: provide detailed run specific information and stop'
111       verbosity=1
112       info_flag=1
113    endif
114    if (trim(inline_options).eq.'-i2') then
115       print*, 'Info mode: provide more detailed run specific information and stop'
116       verbosity=2
117       info_flag=1
118    endif
119  endif
120           
121  if (verbosity.gt.0) then
122    write(*,*) 'call readpaths'
123  endif
124  call readpaths(pathfile)
125 
126  if (verbosity.gt.1) then !show clock info
127     !print*,'length(4)',length(4)
128     !count=0,count_rate=1000
129     CALL SYSTEM_CLOCK(count_clock0, count_rate, count_max)
130     !WRITE(*,*) 'SYSTEM_CLOCK',count, count_rate, count_max
131     !WRITE(*,*) 'SYSTEM_CLOCK, count_clock0', count_clock0
132     !WRITE(*,*) 'SYSTEM_CLOCK, count_rate', count_rate
133     !WRITE(*,*) 'SYSTEM_CLOCK, count_max', count_max
134  endif
135
136  ! Read the user specifications for the current model run
137  !*******************************************************
138
139  if (verbosity.gt.0) then
140    write(*,*) 'call readcommand'
141  endif
142  call readcommand
143  if (verbosity.gt.0) then
144    write(*,*) '    ldirect=', ldirect
145    write(*,*) '    ibdate,ibtime=',ibdate,ibtime
146    write(*,*) '    iedate,ietime=', iedate,ietime
147    if (verbosity.gt.1) then   
148      CALL SYSTEM_CLOCK(count_clock, count_rate, count_max)
149      write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max
150    endif     
151  endif
152
153  ! Read the age classes to be used
154  !********************************
155  if (verbosity.gt.0) then
156    write(*,*) 'call readageclasses'
157  endif
158  call readageclasses
159
160  if (verbosity.gt.1) then   
161    CALL SYSTEM_CLOCK(count_clock, count_rate, count_max)
162    write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max
163  endif     
164
165  ! Read, which wind fields are available within the modelling period
166  !******************************************************************
167
168  if (verbosity.gt.0) then
169    write(*,*) 'call readavailable'
170  endif 
171  call readavailable
172
173  ! If nested wind fields are used, allocate arrays
174  !************************************************
175
176  if (verbosity.gt.0) then
177    write(*,*) 'call com_mod_allocate_nests'
178  endif
179  call com_mod_allocate_nests
180
181  ! Read the model grid specifications,
182  ! both for the mother domain and eventual nests
183  !**********************************************
184 
185  if (verbosity.gt.0) then
186     write(*,*) 'call gridcheck'
187  endif
188
189  call gridcheck
190
191  if (verbosity.gt.1) then   
192    CALL SYSTEM_CLOCK(count_clock, count_rate, count_max)
193    write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max
194  endif     
195
196  if (verbosity.gt.0) then
197    write(*,*) 'call gridcheck_nests'
198  endif 
199  call gridcheck_nests
200
201  ! Read the output grid specifications
202  !************************************
203
204  if (verbosity.gt.0) then
205    write(*,*) 'call readoutgrid'
206  endif
207
208  call readoutgrid
209
210  if (nested_output.eq.1) then
211    call readoutgrid_nest
212    if (verbosity.gt.0) then
213      write(*,*) '# readoutgrid_nest'
214    endif
215  endif
216
217  ! Read the receptor points for which extra concentrations are to be calculated
218  !*****************************************************************************
219
220  if (verbosity.eq.1) then
221     print*,'call readreceptors'
222  endif
223  call readreceptors
224
225  ! Read the physico-chemical species property table
226  !*************************************************
227  !SEC: now only needed SPECIES are read in readreleases.f
228  !call readspecies
229
230
231  ! Read the landuse inventory
232  !***************************
233
234  if (verbosity.gt.0) then
235    print*,'call readlanduse'
236  endif
237  call readlanduse
238
239  ! Assign fractional cover of landuse classes to each ECMWF grid point
240  !********************************************************************
241
242  if (verbosity.gt.0) then
243    print*,'call assignland'
244  endif
245  call assignland
246
247  ! Read the coordinates of the release locations
248  !**********************************************
249
250  if (verbosity.gt.0) then
251    print*,'call readreleases'
252  endif
253  call readreleases
254
255  ! Read and compute surface resistances to dry deposition of gases
256  !****************************************************************
257
258  if (verbosity.gt.0) then
259    print*,'call readdepo'
260  endif
261  call readdepo
262
263  ! Convert the release point coordinates from geografical to grid coordinates
264  !***************************************************************************
265
266  call coordtrafo 
267  if (verbosity.gt.0) then
268    print*,'call coordtrafo'
269  endif
270
271  ! Initialize all particles to non-existent
272  !*****************************************
273
274  if (verbosity.gt.0) then
275    print*,'Initialize all particles to non-existent'
276  endif
277  do j=1,maxpart
278    itra1(j)=-999999999
279  end do
280
281  ! For continuation of previous run, read in particle positions
282  !*************************************************************
283
284  if (ipin.eq.1) then
285    if (verbosity.gt.0) then
286      print*,'call readpartpositions'
287    endif
288    call readpartpositions
289  else
290    if (verbosity.gt.0) then
291      print*,'numpart=0, numparticlecount=0'
292    endif   
293    numpart=0
294    numparticlecount=0
295  endif
296
297  ! Calculate volume, surface area, etc., of all output grid cells
298  ! Allocate fluxes and OHfield if necessary
299  !***************************************************************
300
301  if (verbosity.gt.0) then
302    print*,'call outgrid_init'
303  endif
304  call outgrid_init
305  if (nested_output.eq.1) call outgrid_init_nest
306
307  ! Read the OH field
308  !******************
309
310  if (OHREA.eqv..TRUE.) then
311    if (verbosity.gt.0) then
312      print*,'call readOHfield'
313    endif
314    call readOHfield
315  endif
316
317  ! Write basic information on the simulation to a file "header"
318  ! and open files that are to be kept open throughout the simulation
319  !******************************************************************
320
321  if (lnetcdfout.eq.1) then
322    call writeheader_netcdf(lnest=.false.)
323  else
324    call writeheader
325  end if
326
327  if (nested_output.eq.1) then
328    if (lnetcdfout.eq.1) then
329      call writeheader_netcdf(lnest=.true.)
330    else
331      call writeheader_nest
332    endif
333  endif
334
335  if (verbosity.gt.0) then
336    print*,'call writeheader'
337  endif
338
339  call writeheader
340  ! FLEXPART 9.2 ticket ?? write header in ASCII format
341  call writeheader_txt
342  !if (nested_output.eq.1) call writeheader_nest
343  if (nested_output.eq.1.and.surf_only.ne.1) call writeheader_nest
344  if (nested_output.eq.1.and.surf_only.eq.1) call writeheader_nest_surf
345  if (nested_output.ne.1.and.surf_only.eq.1) call writeheader_surf
346
347  !open(unitdates,file=path(2)(1:length(2))//'dates')
348
349  if (verbosity.gt.0) then
350    print*,'call openreceptors'
351  endif
352  call openreceptors
353  if ((iout.eq.4).or.(iout.eq.5)) call openouttraj
354
355  ! Releases can only start and end at discrete times (multiples of lsynctime)
356  !***************************************************************************
357
358  if (verbosity.gt.0) then
359    print*,'discretize release times'
360  endif
361  do i=1,numpoint
362    ireleasestart(i)=nint(real(ireleasestart(i))/real(lsynctime))*lsynctime
363    ireleaseend(i)=nint(real(ireleaseend(i))/real(lsynctime))*lsynctime
364  end do
365
366  ! Initialize cloud-base mass fluxes for the convection scheme
367  !************************************************************
368
369  if (verbosity.gt.0) then
370    print*,'Initialize cloud-base mass fluxes for the convection scheme'
371  endif
372
373  do jy=0,nymin1
374    do ix=0,nxmin1
375      cbaseflux(ix,jy)=0.
376    end do
377  end do
378  do inest=1,numbnests
379    do jy=0,nyn(inest)-1
380      do ix=0,nxn(inest)-1
381        cbasefluxn(ix,jy,inest)=0.
382      end do
383    end do
384  end do
385
386  ! Calculate particle trajectories
387  !********************************
388
389  if (verbosity.gt.0) then
390     if (verbosity.gt.1) then   
391       CALL SYSTEM_CLOCK(count_clock, count_rate, count_max)
392       write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max
393     endif
394     if (info_flag.eq.1) then
395       print*, 'info only mode (stop)'   
396       stop
397     endif
398     print*,'call timemanager'
399  endif
400
401  call timemanager
402
403! NIK 16.02.2005
404  write(*,*) '**********************************************'
405  write(*,*) 'Total number of occurences of below-cloud scavenging', tot_blc_count
406  write(*,*) 'Total number of occurences of in-cloud    scavenging', tot_inc_count
407  write(*,*) '**********************************************'
408
409  write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
410       &XPART MODEL RUN!'
411
412end program flexpart
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG