source: flexpart.git/src/FLEXPART_MPI.f90 @ fddc6ec

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

Added a version number string for flexpart version, that will also be used for SPECIES file names. TODO: update unit tests, species and readspecies (master branch?)

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