source: flexpart.git/src/readpartpositions.f90 @ 02095e3

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

Updated checks and warning messages for wet deposition parameters given in SPECIES files

  • Property mode set to 100644
File size: 6.3 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
22subroutine readpartpositions
23
24  !*****************************************************************************
25  !                                                                            *
26  !   This routine opens the particle dump file and reads all the particle     *
27  !   positions from a previous run to initialize the current run.             *
28  !                                                                            *
29  !                                                                            *
30  !     Author: A. Stohl                                                       *
31  !                                                                            *
32  !     24 March 2000                                                          *
33  !                                                                            *
34  !*****************************************************************************
35  !                                                                            *
36  ! Variables:                                                                 *
37  !                                                                            *
38  !*****************************************************************************
39
40  use par_mod
41  use com_mod
42  use random_mod
43
44  implicit none
45
46  integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix
47  integer :: id1,id2,it1,it2
48  real :: xlonin,ylatin,topo,hmixi,pvi,qvi,rhoi,tri,tti
49  character :: specin*7
50  real(kind=dp) :: julin,julpartin,juldate
51
52  integer :: idummy = -8
53
54  numparticlecount=0
55
56  ! Open header file of dumped particle data
57  !*****************************************
58
59  open(unitpartin,file=path(2)(1:length(2))//'header', &
60       form='unformatted',err=998)
61
62  read(unitpartin) ibdatein,ibtimein
63  read(unitpartin)
64  read(unitpartin)
65
66  read(unitpartin)
67  read(unitpartin)
68  read(unitpartin) nspecin
69  nspecin=nspecin/3
70  if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997
71
72  do i=1,nspecin
73    read(unitpartin)
74    read(unitpartin)
75    read(unitpartin) j,specin
76    if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996
77  end do
78
79  read(unitpartin) numpointin
80  if (numpointin.ne.numpoint) goto 995
81999 continue
82  do i=1,numpointin
83    read(unitpartin)
84    read(unitpartin)
85    read(unitpartin)
86    read(unitpartin)
87    do j=1,nspec
88      read(unitpartin)
89      read(unitpartin)
90      read(unitpartin)
91    end do
92  end do
93  read(unitpartin)
94  read(unitpartin)
95
96  do ix=0,numxgrid-1
97    read(unitpartin)
98  end do
99
100
101  ! Open data file of dumped particle data
102  !***************************************
103
104  close(unitpartin)
105  open(unitpartin,file=path(2)(1:length(2))//'partposit_end', &
106       form='unformatted',err=998)
107 
108
109100 read(unitpartin,end=99) itimein
110  i=0
111200 i=i+1
112  read(unitpartin) npoint(i),xlonin,ylatin,ztra1(i),itramem(i), &
113       topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1(i,j),j=1,nspec)
114 
115  if (xlonin.eq.-9999.9) goto 100
116  xtra1(i)=(xlonin-xlon0)/dx
117  ytra1(i)=(ylatin-ylat0)/dy
118  numparticlecount=max(numparticlecount,npoint(i))
119  goto 200
120
12199 numpart=i-1
122
123  close(unitpartin)
124
125  julin=juldate(ibdatein,ibtimein)+real(itimein,kind=dp)/86400._dp
126  if (abs(julin-bdate).gt.1.e-5) goto 994
127  do i=1,numpart
128    julpartin=juldate(ibdatein,ibtimein)+ &
129         real(itramem(i),kind=dp)/86400._dp
130    nclass(i)=min(int(ran1(idummy)*real(nclassunc))+1, &
131         nclassunc)
132    idt(i)=mintime
133    itra1(i)=0
134    itramem(i)=nint((julpartin-bdate)*86400.)
135    itrasplit(i)=ldirect*itsplit
136  end do
137
138  return
139
140
141994   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
142  write(*,*) ' #### ENDING TIME OF PREVIOUS MODEL RUN DOES   #### '
143  write(*,*) ' #### NOT AGREE WITH STARTING TIME OF THIS RUN.#### '
144  call caldate(julin,id1,it1)
145  call caldate(bdate,id2,it2)
146  write(*,*) 'julin: ',julin,id1,it1
147  write(*,*) 'bdate: ',bdate,id2,it2
148  stop
149
150!995   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
151995   write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
152  write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT     #### '
153  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
154!  stop
155  goto 999
156
157996   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
158  write(*,*) ' #### SPECIES NAMES TO BE READ IN DO NOT       #### '
159  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
160  stop
161
162997   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
163  write(*,*) ' #### THE NUMBER OF SPECIES TO BE READ IN DOES #### '
164  write(*,*) ' #### NOT AGREE WITH CURRENT SETTINGS!         #### '
165  stop
166
167998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
168  write(*,*) ' #### '//path(2)(1:length(2))//'grid'//' #### '
169  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
170  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
171  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
172  stop
173
174end subroutine readpartpositions
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG