source: flexpart.git/src/readpartpositions.f90

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

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 5.0 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine readpartpositions
5
6  !*****************************************************************************
7  !                                                                            *
8  !   This routine opens the particle dump file and reads all the particle     *
9  !   positions from a previous run to initialize the current run.             *
10  !                                                                            *
11  !                                                                            *
12  !     Author: A. Stohl                                                       *
13  !                                                                            *
14  !     24 March 2000                                                          *
15  !                                                                            *
16  !*****************************************************************************
17  !                                                                            *
18  ! Variables:                                                                 *
19  !                                                                            *
20  !*****************************************************************************
21
22  use par_mod
23  use com_mod
24  use random_mod
25
26  implicit none
27
28  integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix
29  integer :: id1,id2,it1,it2
30  real :: xlonin,ylatin,topo,hmixi,pvi,qvi,rhoi,tri,tti
31  character :: specin*7
32  real(kind=dp) :: julin,julpartin,juldate
33
34  integer :: idummy = -8
35
36  numparticlecount=0
37
38  ! Open header file of dumped particle data
39  !*****************************************
40
41  open(unitpartin,file=path(2)(1:length(2))//'header', &
42       form='unformatted',err=998)
43
44  read(unitpartin) ibdatein,ibtimein
45  read(unitpartin)
46  read(unitpartin)
47
48  read(unitpartin)
49  read(unitpartin)
50  read(unitpartin) nspecin
51  nspecin=nspecin/3
52  if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997
53
54  do i=1,nspecin
55    read(unitpartin)
56    read(unitpartin)
57    read(unitpartin) j,specin
58    if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996
59  end do
60
61  read(unitpartin) numpointin
62  if (numpointin.ne.numpoint) goto 995
63999 continue
64  do i=1,numpointin
65    read(unitpartin)
66    read(unitpartin)
67    read(unitpartin)
68    read(unitpartin)
69    do j=1,nspec
70      read(unitpartin)
71      read(unitpartin)
72      read(unitpartin)
73    end do
74  end do
75  read(unitpartin)
76  read(unitpartin)
77
78  do ix=0,numxgrid-1
79    read(unitpartin)
80  end do
81
82
83  ! Open data file of dumped particle data
84  !***************************************
85
86  close(unitpartin)
87  open(unitpartin,file=path(2)(1:length(2))//'partposit_end', &
88       form='unformatted',err=998)
89 
90
91100 read(unitpartin,end=99) itimein
92  i=0
93200 i=i+1
94  read(unitpartin) npoint(i),xlonin,ylatin,ztra1(i),itramem(i), &
95       topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1(i,j),j=1,nspec)
96 
97  if (xlonin.eq.-9999.9) goto 100
98  xtra1(i)=(xlonin-xlon0)/dx
99  ytra1(i)=(ylatin-ylat0)/dy
100  numparticlecount=max(numparticlecount,npoint(i))
101  goto 200
102
10399 numpart=i-1
104
105  close(unitpartin)
106
107  julin=juldate(ibdatein,ibtimein)+real(itimein,kind=dp)/86400._dp
108  if (abs(julin-bdate).gt.1.e-5) goto 994
109  do i=1,numpart
110    julpartin=juldate(ibdatein,ibtimein)+ &
111         real(itramem(i),kind=dp)/86400._dp
112    nclass(i)=min(int(ran1(idummy)*real(nclassunc))+1, &
113         nclassunc)
114    idt(i)=mintime
115    itra1(i)=0
116    itramem(i)=nint((julpartin-bdate)*86400.)
117    itrasplit(i)=ldirect*itsplit
118  end do
119
120  return
121
122
123994   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
124  write(*,*) ' #### ENDING TIME OF PREVIOUS MODEL RUN DOES   #### '
125  write(*,*) ' #### NOT AGREE WITH STARTING TIME OF THIS RUN.#### '
126  call caldate(julin,id1,it1)
127  call caldate(bdate,id2,it2)
128  write(*,*) 'julin: ',julin,id1,it1
129  write(*,*) 'bdate: ',bdate,id2,it2
130  stop
131
132!995   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
133995   write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
134  write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT     #### '
135  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
136!  stop
137  goto 999
138
139996   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
140  write(*,*) ' #### SPECIES NAMES TO BE READ IN DO NOT       #### '
141  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
142  stop
143
144997   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
145  write(*,*) ' #### THE NUMBER OF SPECIES TO BE READ IN DOES #### '
146  write(*,*) ' #### NOT AGREE WITH CURRENT SETTINGS!         #### '
147  stop
148
149998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
150  write(*,*) ' #### '//path(2)(1:length(2))//'grid'//' #### '
151  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
152  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
153  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
154  stop
155
156end subroutine readpartpositions
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG