[e200b7a] | 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 | |
---|
| 22 | subroutine 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 | |
---|
| 43 | implicit none |
---|
| 44 | |
---|
| 45 | integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix |
---|
| 46 | integer :: id1,id2,it1,it2 |
---|
| 47 | real :: xlonin,ylatin,ran1,topo,hmixi,pvi,qvi,rhoi,tri,tti |
---|
| 48 | character :: specin*7 |
---|
| 49 | real(kind=dp) :: julin,julpartin,juldate |
---|
| 50 | |
---|
| 51 | integer :: idummy = -8 |
---|
| 52 | |
---|
| 53 | numparticlecount=0 |
---|
| 54 | |
---|
| 55 | ! Open header file of dumped particle data |
---|
| 56 | !***************************************** |
---|
| 57 | |
---|
| 58 | open(unitpartin,file=path(2)(1:length(2))//'header', & |
---|
| 59 | form='unformatted',err=998) |
---|
| 60 | |
---|
| 61 | read(unitpartin) ibdatein,ibtimein |
---|
| 62 | read(unitpartin) |
---|
| 63 | read(unitpartin) |
---|
| 64 | |
---|
| 65 | read(unitpartin) |
---|
| 66 | read(unitpartin) |
---|
| 67 | read(unitpartin) nspecin |
---|
| 68 | nspecin=nspecin/3 |
---|
| 69 | if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997 |
---|
| 70 | |
---|
| 71 | do i=1,nspecin |
---|
| 72 | read(unitpartin) |
---|
| 73 | read(unitpartin) |
---|
| 74 | read(unitpartin) j,specin |
---|
| 75 | if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996 |
---|
| 76 | end do |
---|
| 77 | |
---|
| 78 | read(unitpartin) numpointin |
---|
[f13406c] | 79 | if (numpointin.ne.numpoint) goto 995 |
---|
| 80 | 999 continue |
---|
[e200b7a] | 81 | do i=1,numpointin |
---|
| 82 | read(unitpartin) |
---|
| 83 | read(unitpartin) |
---|
| 84 | read(unitpartin) |
---|
| 85 | read(unitpartin) |
---|
| 86 | do j=1,nspec |
---|
| 87 | read(unitpartin) |
---|
| 88 | read(unitpartin) |
---|
| 89 | read(unitpartin) |
---|
| 90 | end do |
---|
| 91 | end do |
---|
| 92 | read(unitpartin) |
---|
| 93 | read(unitpartin) |
---|
| 94 | |
---|
| 95 | do ix=0,numxgrid-1 |
---|
| 96 | read(unitpartin) |
---|
| 97 | end do |
---|
| 98 | |
---|
| 99 | |
---|
| 100 | ! Open data file of dumped particle data |
---|
| 101 | !*************************************** |
---|
| 102 | |
---|
| 103 | close(unitpartin) |
---|
| 104 | open(unitpartin,file=path(2)(1:length(2))//'partposit_end', & |
---|
| 105 | form='unformatted',err=998) |
---|
| 106 | |
---|
| 107 | |
---|
| 108 | 100 read(unitpartin,end=99) itimein |
---|
| 109 | i=0 |
---|
| 110 | 200 i=i+1 |
---|
| 111 | read(unitpartin) npoint(i),xlonin,ylatin,ztra1(i),itramem(i), & |
---|
| 112 | topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1(i,j),j=1,nspec) |
---|
| 113 | |
---|
| 114 | if (xlonin.eq.-9999.9) goto 100 |
---|
| 115 | xtra1(i)=(xlonin-xlon0)/dx |
---|
| 116 | ytra1(i)=(ylatin-ylat0)/dy |
---|
| 117 | numparticlecount=max(numparticlecount,npoint(i)) |
---|
| 118 | goto 200 |
---|
| 119 | |
---|
| 120 | 99 numpart=i-1 |
---|
| 121 | |
---|
| 122 | close(unitpartin) |
---|
| 123 | |
---|
| 124 | julin=juldate(ibdatein,ibtimein)+real(itimein,kind=dp)/86400._dp |
---|
| 125 | if (abs(julin-bdate).gt.1.e-5) goto 994 |
---|
| 126 | do i=1,numpart |
---|
| 127 | julpartin=juldate(ibdatein,ibtimein)+ & |
---|
| 128 | real(itramem(i),kind=dp)/86400._dp |
---|
| 129 | nclass(i)=min(int(ran1(idummy)*real(nclassunc))+1, & |
---|
| 130 | nclassunc) |
---|
| 131 | idt(i)=mintime |
---|
| 132 | itra1(i)=0 |
---|
| 133 | itramem(i)=nint((julpartin-bdate)*86400.) |
---|
| 134 | itrasplit(i)=ldirect*itsplit |
---|
| 135 | end do |
---|
| 136 | |
---|
| 137 | return |
---|
| 138 | |
---|
| 139 | |
---|
| 140 | 994 write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### ' |
---|
| 141 | write(*,*) ' #### ENDING TIME OF PREVIOUS MODEL RUN DOES #### ' |
---|
| 142 | write(*,*) ' #### NOT AGREE WITH STARTING TIME OF THIS RUN.#### ' |
---|
| 143 | call caldate(julin,id1,it1) |
---|
| 144 | call caldate(bdate,id2,it2) |
---|
| 145 | write(*,*) 'julin: ',julin,id1,it1 |
---|
| 146 | write(*,*) 'bdate: ',bdate,id2,it2 |
---|
| 147 | stop |
---|
| 148 | |
---|
[f13406c] | 149 | !995 write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### ' |
---|
| 150 | 995 write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### ' |
---|
| 151 | write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT #### ' |
---|
| 152 | write(*,*) ' #### AGREE WITH CURRENT SETTINGS! #### ' |
---|
| 153 | ! stop |
---|
| 154 | goto 999 |
---|
[e200b7a] | 155 | |
---|
| 156 | 996 write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### ' |
---|
| 157 | write(*,*) ' #### SPECIES NAMES TO BE READ IN DO NOT #### ' |
---|
| 158 | write(*,*) ' #### AGREE WITH CURRENT SETTINGS! #### ' |
---|
| 159 | stop |
---|
| 160 | |
---|
| 161 | 997 write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### ' |
---|
| 162 | write(*,*) ' #### THE NUMBER OF SPECIES TO BE READ IN DOES #### ' |
---|
| 163 | write(*,*) ' #### NOT AGREE WITH CURRENT SETTINGS! #### ' |
---|
| 164 | stop |
---|
| 165 | |
---|
| 166 | 998 write(*,*) ' #### FLEXPART MODEL ERROR! THE FILE #### ' |
---|
| 167 | write(*,*) ' #### '//path(2)(1:length(2))//'grid'//' #### ' |
---|
| 168 | write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS #### ' |
---|
| 169 | write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### ' |
---|
| 170 | write(*,*) ' #### THE PROGRAM AGAIN. #### ' |
---|
| 171 | stop |
---|
| 172 | |
---|
| 173 | end subroutine readpartpositions |
---|