source: flexpart.git/src/partoutput_mpi.f90 @ 3481cc1

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

move license from headers to a different file

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