source: trunk/src/partoutput.f90

Last change on this file was 4, checked in by mlanger, 11 years ago
File size: 7.1 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 partoutput(itime)
23  !                        i
24  !*****************************************************************************
25  !                                                                            *
26  !     Dump all particle positions                                            *
27  !                                                                            *
28  !     Author: A. Stohl                                                       *
29  !                                                                            *
30  !     12 March 1999                                                          *
31  !                                                                            *
32  !*****************************************************************************
33  !                                                                            *
34  ! Variables:                                                                 *
35  !                                                                            *
36  !*****************************************************************************
37
38  use par_mod
39  use com_mod
40
41  implicit none
42
43  real(kind=dp) :: jul
44  integer :: itime,i,j,jjjjmmdd,ihmmss
45  integer :: ix,jy,ixp,jyp,indexh,m,il,ind,indz,indzp
46  real :: xlon,ylat
47  real :: dt1,dt2,dtt,ddx,ddy,rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
48  real :: topo,hm(2),hmixi,pv1(2),pvprof(2),pvi,qv1(2),qvprof(2),qvi
49  real :: tt1(2),ttprof(2),tti,rho1(2),rhoprof(2),rhoi
50  real :: tr(2),tri
51  character :: adate*8,atime*6
52
53
54  ! Determine current calendar date, needed for the file name
55  !**********************************************************
56
57  jul=bdate+real(itime,kind=dp)/86400._dp
58  call caldate(jul,jjjjmmdd,ihmmss)
59  write(adate,'(i8.8)') jjjjmmdd
60  write(atime,'(i6.6)') ihmmss
61
62
63  ! Some variables needed for temporal interpolation
64  !*************************************************
65
66  dt1=real(itime-memtime(1))
67  dt2=real(memtime(2)-itime)
68  dtt=1./(dt1+dt2)
69
70  ! Open output file and write the output
71  !**************************************
72
73  if (ipout.eq.1) then
74    open(unitpartout,file=path(2)(1:length(2))//'partposit_'//adate// &
75         atime,form='unformatted')
76  else
77    open(unitpartout,file=path(2)(1:length(2))//'partposit_end', &
78         form='unformatted')
79  endif
80
81  ! Write current time to file
82  !***************************
83
84  write(unitpartout) itime
85  do i=1,numpart
86
87  ! Take only valid particles
88  !**************************
89
90    if (itra1(i).eq.itime) then
91      xlon=xlon0+xtra1(i)*dx
92      ylat=ylat0+ytra1(i)*dy
93
94  !*****************************************************************************
95  ! Interpolate several variables (PV, specific humidity, etc.) to particle position
96  !*****************************************************************************
97
98      ix=xtra1(i)
99      jy=ytra1(i)
100      ixp=ix+1
101      jyp=jy+1
102      ddx=xtra1(i)-real(ix)
103      ddy=ytra1(i)-real(jy)
104      rddx=1.-ddx
105      rddy=1.-ddy
106      p1=rddx*rddy
107      p2=ddx*rddy
108      p3=rddx*ddy
109      p4=ddx*ddy
110
111  ! Topography
112  !***********
113
114      topo=p1*oro(ix ,jy) &
115           + p2*oro(ixp,jy) &
116           + p3*oro(ix ,jyp) &
117           + p4*oro(ixp,jyp)
118
119  ! Potential vorticity, specific humidity, temperature, and density
120  !*****************************************************************
121
122      do il=2,nz
123        if (height(il).gt.ztra1(i)) then
124          indz=il-1
125          indzp=il
126          goto 6
127        endif
128      end do
1296     continue
130
131      dz1=ztra1(i)-height(indz)
132      dz2=height(indzp)-ztra1(i)
133      dz=1./(dz1+dz2)
134
135
136      do ind=indz,indzp
137        do m=1,2
138          indexh=memind(m)
139
140  ! Potential vorticity
141          pv1(m)=p1*pv(ix ,jy ,ind,indexh) &
142               +p2*pv(ixp,jy ,ind,indexh) &
143               +p3*pv(ix ,jyp,ind,indexh) &
144               +p4*pv(ixp,jyp,ind,indexh)
145  ! Specific humidity
146          qv1(m)=p1*qv(ix ,jy ,ind,indexh) &
147               +p2*qv(ixp,jy ,ind,indexh) &
148               +p3*qv(ix ,jyp,ind,indexh) &
149               +p4*qv(ixp,jyp,ind,indexh)
150  ! Temperature
151          tt1(m)=p1*tt(ix ,jy ,ind,indexh) &
152               +p2*tt(ixp,jy ,ind,indexh) &
153               +p3*tt(ix ,jyp,ind,indexh) &
154               +p4*tt(ixp,jyp,ind,indexh)
155  ! Density
156          rho1(m)=p1*rho(ix ,jy ,ind,indexh) &
157               +p2*rho(ixp,jy ,ind,indexh) &
158               +p3*rho(ix ,jyp,ind,indexh) &
159               +p4*rho(ixp,jyp,ind,indexh)
160        end do
161        pvprof(ind-indz+1)=(pv1(1)*dt2+pv1(2)*dt1)*dtt
162        qvprof(ind-indz+1)=(qv1(1)*dt2+qv1(2)*dt1)*dtt
163        ttprof(ind-indz+1)=(tt1(1)*dt2+tt1(2)*dt1)*dtt
164        rhoprof(ind-indz+1)=(rho1(1)*dt2+rho1(2)*dt1)*dtt
165      end do
166      pvi=(dz1*pvprof(2)+dz2*pvprof(1))*dz
167      qvi=(dz1*qvprof(2)+dz2*qvprof(1))*dz
168      tti=(dz1*ttprof(2)+dz2*ttprof(1))*dz
169      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
170
171  ! Tropopause and PBL height
172  !**************************
173
174      do m=1,2
175        indexh=memind(m)
176
177  ! Tropopause
178        tr(m)=p1*tropopause(ix ,jy ,1,indexh) &
179             + p2*tropopause(ixp,jy ,1,indexh) &
180             + p3*tropopause(ix ,jyp,1,indexh) &
181             + p4*tropopause(ixp,jyp,1,indexh)
182
183  ! PBL height
184        hm(m)=p1*hmix(ix ,jy ,1,indexh) &
185             + p2*hmix(ixp,jy ,1,indexh) &
186             + p3*hmix(ix ,jyp,1,indexh) &
187             + p4*hmix(ixp,jyp,1,indexh)
188      end do
189
190      hmixi=(hm(1)*dt2+hm(2)*dt1)*dtt
191      tri=(tr(1)*dt2+tr(2)*dt1)*dtt
192
193
194  ! Write the output
195  !*****************
196
197      write(unitpartout) npoint(i),xlon,ylat,ztra1(i), &
198           itramem(i),topo,pvi,qvi,rhoi,hmixi,tri,tti, &
199           (xmass1(i,j),j=1,nspec)
200    endif
201  end do
202  write(unitpartout) -99999,-9999.9,-9999.9,-9999.9,-99999, &
203       -9999.9,-9999.9,-9999.9,-9999.9,-9999.9,-9999.9,-9999.9, &
204       (-9999.9,j=1,nspec)
205
206
207  close(unitpartout)
208
209end subroutine partoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG