partoutput_average_mpi.f90 7.11 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
subroutine partoutput_average(itime,irec)
  !                             i
  !*****************************************************************************
  !                                                                            *
  !     Dump all particle positions                                            *
  !                                                                            *
  !     Author: A. Stohl                                                       *
  !                                                                            *
  !     12 March 1999                                                          *
  !     03/2019 AST: Version for MPI                                           *
  !                  processes sequentially access and append data to file     *
  !                                                                            *
  !                                                                            *
  !*****************************************************************************
  !                                                                            *
  ! Variables:                                                                 *
  !                                                                            *
  !*****************************************************************************

  use par_mod
  use com_mod
  use mpi_mod


  implicit none

  real(kind=dp) :: jul
  integer :: itime,i,j,jjjjmmdd,ihmmss,irec
  integer :: ix,jy,ixp,jyp,indexh,m,il,ind,indz,indzp
  real :: xlon,ylat
  real :: dt1,dt2,dtt,ddx,ddy,rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
  real :: topo,hm(2),hmixi,pv1(2),pvprof(2),pvi,qv1(2),qvprof(2),qvi
  real :: tt1(2),ttprof(2),tti,rho1(2),rhoprof(2),rhoi
  real :: tr(2),tri,zlim
  character :: adate*8,atime*6
  character(LEN=8) :: file_stat='OLD'

  integer(kind=2) :: ishort_xlon(maxpart_mpi),ishort_ylat(maxpart_mpi),ishort_z(maxpart_mpi)
  integer(kind=2) :: ishort_topo(maxpart_mpi),ishort_tro(maxpart_mpi),ishort_hmix(maxpart_mpi)
  integer(kind=2) :: ishort_pv(maxpart_mpi),ishort_rho(maxpart_mpi),ishort_qv(maxpart_mpi)
  integer(kind=2) :: ishort_tt(maxpart_mpi),ishort_uu(maxpart_mpi),ishort_vv(maxpart_mpi)
  integer(kind=2) :: ishort_energy(maxpart_mpi)

  ! MPI root process creates/overwrites the file, other processes append to it
  if (lroot) file_stat='REPLACE'

  ! Determine current calendar date, needed for the file name
  !**********************************************************

  jul=bdate+real(itime,kind=dp)/86400._dp
  call caldate(jul,jjjjmmdd,ihmmss)
  write(adate,'(i8.8)') jjjjmmdd
  write(atime,'(i6.6)') ihmmss


  ! Some variables needed for temporal interpolation
  !*************************************************

  dt1=real(itime-memtime(1))
  dt2=real(memtime(2)-itime)
  dtt=1./(dt1+dt2)

  ! Open output file and write the output
  !**************************************

  open(unitpartout_average,file=path(2)(1:length(2))//'partposit_average_'//adate// &
       atime,form='unformatted',access='direct',status=file_stat,recl=24)


  ! Write current time to file
  !***************************

!  write(unitpartout_average) itime,numpart
  do i=1,numpart

  ! Take only valid particles
  !**************************

    if (itra1(i).eq.itime) then
      part_av_topo(i)=part_av_topo(i)/float(npart_av(i))
      part_av_z(i)=part_av_z(i)/float(npart_av(i))
      part_av_pv(i)=part_av_pv(i)/float(npart_av(i))
      part_av_qv(i)=part_av_qv(i)/float(npart_av(i))
      part_av_tt(i)=part_av_tt(i)/float(npart_av(i))
      part_av_uu(i)=part_av_uu(i)/float(npart_av(i))
      part_av_vv(i)=part_av_vv(i)/float(npart_av(i))
      part_av_rho(i)=part_av_rho(i)/float(npart_av(i))
      part_av_tro(i)=part_av_tro(i)/float(npart_av(i))
      part_av_hmix(i)=part_av_hmix(i)/float(npart_av(i))
      part_av_energy(i)=part_av_energy(i)/float(npart_av(i))

      part_av_cartx(i)=part_av_cartx(i)/float(npart_av(i))
      part_av_carty(i)=part_av_carty(i)/float(npart_av(i))
      part_av_cartz(i)=part_av_cartz(i)/float(npart_av(i))

! Project Cartesian coordinates back onto Earth's surface
! #######################################################
      xlon=atan2(part_av_cartx(i),-1.*part_av_carty(i))
      ylat=atan2(part_av_cartz(i),sqrt(part_av_cartx(i)*part_av_cartx(i)+ &
           part_av_carty(i)*part_av_carty(i)))
      xlon=xlon/pi180
      ylat=ylat/pi180


  ! Convert all data to integer*2 variables (from -32768 to 32767) for compressed output
  !*************************************************************************************

      if (xlon.gt.180.) xlon=xlon-360.
      if (xlon.lt.-180.) xlon=xlon+360.
      ishort_xlon(i)=nint(xlon*180.)
      ishort_ylat(i)=nint(ylat*360.)

      zlim=(part_av_z(i)*2.-32000.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_z(i)=nint(zlim)

      zlim=(part_av_topo(i)*2.-32000.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_topo(i)=nint(zlim)

      zlim=(part_av_tro(i)*2.-32000.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_tro(i)=nint(zlim)

      zlim=(part_av_hmix(i)*2.-32000.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_hmix(i)=nint(zlim)

      zlim=(part_av_rho(i)*20000.-32000.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_rho(i)=nint(zlim)

      zlim=(part_av_qv(i)*1000000.-32000.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_qv(i)=nint(zlim)

      zlim=(part_av_pv(i)*100.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_pv(i)=nint(zlim)

      zlim=((part_av_tt(i)-273.15))*300.
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_tt(i)=nint(zlim)

      zlim=(part_av_uu(i)*200.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_uu(i)=nint(zlim)

      zlim=(part_av_vv(i)*200.)
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_vv(i)=nint(zlim)

      zlim=(part_av_energy(i)-300000.)/30.
      zlim=min(zlim,32766.)
      zlim=max(zlim,-32766.)
      ishort_energy(i)=nint(zlim)


  ! Turn on for readable test output
  !*********************************

!        write(119,*) itime,i,xlon,ylat,part_av_z(i),part_av_topo(i),part_av_tro(i), &
!        part_av_hmix(i),part_av_rho(i),part_av_qv(i),part_av_pv(i),part_av_tt(i), &
!        ishort_uu(i),ishort_vv(i)
    endif

! Re-initialize averages, and set number of steps to zero
    npart_av(i)=0
    part_av_topo(i)=0.
    part_av_z(i)=0.
    part_av_cartx(i)=0.
    part_av_carty(i)=0.
    part_av_cartz(i)=0.
    part_av_pv(i)=0.
    part_av_qv(i)=0.
    part_av_tt(i)=0.
    part_av_uu(i)=0.
    part_av_vv(i)=0.
    part_av_rho(i)=0.
    part_av_tro(i)=0.
    part_av_hmix(i)=0.
    part_av_energy(i)=0.
    
  end do

  ! Write the output
  !*****************

  do i=1,numpart
    if (itra1(i).eq.itime) then
      write(unitpartout_average,rec=irec+i) ishort_xlon(i),ishort_ylat(i),ishort_z(i), &
           ishort_topo(i),ishort_tro(i),ishort_hmix(i),ishort_rho(i),ishort_qv(i),ishort_pv(i), &
           ishort_tt(i),ishort_uu(i),ishort_vv(i)  ! ,ishort_energy(i)
    endif
  enddo

  close(unitpartout_average)

end subroutine partoutput_average