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