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