source: flexpart.git/src/partoutput_mpi.f90 @ 77778f8

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 77778f8 was d2a5a83, checked in by Espen Sollum ATMOS <eso@…>, 6 years ago

Fixed an issue with allocation in previous commit

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