source: flexpart.git/src/FLEXPART.f90 @ 2753a5c

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 2753a5c was 2753a5c, checked in by Ignacio Pisso <ip@…>, 5 years ago

resolved merging conflicts with GFS branch

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