source: flexpart.git/src/FLEXPART.f90 @ 414a5e5

10.4.1_peseiFPv9.3.1FPv9.3.1b_testingFPv9.3.2GFS_025bugfixes+enhancementsdevfp9.3.1-20161214-nc4grib2nc4_repairrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 414a5e5 was 414a5e5, checked in by Ignacio Pisso <Ignacio.Pisso@…>, 9 years ago

add verbose message to FLEXPART.f90

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