[0c8c7f2] | 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 | |
---|
| 22 | subroutine partoutput_average(itime,irec) |
---|
| 23 | ! i |
---|
| 24 | !***************************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! Dump all particle positions * |
---|
| 27 | ! * |
---|
| 28 | ! Author: A. Stohl * |
---|
| 29 | ! * |
---|
| 30 | ! 12 March 1999 * |
---|
| 31 | ! 03/2019 AST: Version for MPI * |
---|
| 32 | ! processes sequentially access and append data to file * |
---|
| 33 | ! * |
---|
| 34 | ! * |
---|
| 35 | !***************************************************************************** |
---|
| 36 | ! * |
---|
| 37 | ! Variables: * |
---|
| 38 | ! * |
---|
| 39 | !***************************************************************************** |
---|
| 40 | |
---|
| 41 | use par_mod |
---|
| 42 | use com_mod |
---|
| 43 | use mpi_mod |
---|
| 44 | |
---|
| 45 | |
---|
| 46 | implicit none |
---|
| 47 | |
---|
| 48 | real(kind=dp) :: jul |
---|
| 49 | integer :: itime,i,j,jjjjmmdd,ihmmss,irec |
---|
| 50 | integer :: ix,jy,ixp,jyp,indexh,m,il,ind,indz,indzp |
---|
| 51 | real :: xlon,ylat |
---|
| 52 | real :: dt1,dt2,dtt,ddx,ddy,rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz |
---|
| 53 | real :: topo,hm(2),hmixi,pv1(2),pvprof(2),pvi,qv1(2),qvprof(2),qvi |
---|
| 54 | real :: tt1(2),ttprof(2),tti,rho1(2),rhoprof(2),rhoi |
---|
| 55 | real :: tr(2),tri,zlim |
---|
| 56 | character :: adate*8,atime*6 |
---|
| 57 | character(LEN=8) :: file_stat='OLD' |
---|
| 58 | |
---|
| 59 | integer(kind=2) :: ishort_xlon(maxpart_mpi),ishort_ylat(maxpart_mpi),ishort_z(maxpart_mpi) |
---|
| 60 | integer(kind=2) :: ishort_topo(maxpart_mpi),ishort_tro(maxpart_mpi),ishort_hmix(maxpart_mpi) |
---|
| 61 | integer(kind=2) :: ishort_pv(maxpart_mpi),ishort_rho(maxpart_mpi),ishort_qv(maxpart_mpi) |
---|
| 62 | integer(kind=2) :: ishort_tt(maxpart_mpi),ishort_uu(maxpart_mpi),ishort_vv(maxpart_mpi) |
---|
| 63 | integer(kind=2) :: ishort_energy(maxpart_mpi) |
---|
| 64 | |
---|
| 65 | ! MPI root process creates/overwrites the file, other processes append to it |
---|
| 66 | if (lroot) file_stat='REPLACE' |
---|
| 67 | |
---|
| 68 | ! Determine current calendar date, needed for the file name |
---|
| 69 | !********************************************************** |
---|
| 70 | |
---|
| 71 | jul=bdate+real(itime,kind=dp)/86400._dp |
---|
| 72 | call caldate(jul,jjjjmmdd,ihmmss) |
---|
| 73 | write(adate,'(i8.8)') jjjjmmdd |
---|
| 74 | write(atime,'(i6.6)') ihmmss |
---|
| 75 | |
---|
| 76 | |
---|
| 77 | ! Some variables needed for temporal interpolation |
---|
| 78 | !************************************************* |
---|
| 79 | |
---|
| 80 | dt1=real(itime-memtime(1)) |
---|
| 81 | dt2=real(memtime(2)-itime) |
---|
| 82 | dtt=1./(dt1+dt2) |
---|
| 83 | |
---|
| 84 | ! Open output file and write the output |
---|
| 85 | !************************************** |
---|
| 86 | |
---|
| 87 | open(unitpartout_average,file=path(2)(1:length(2))//'partposit_average_'//adate// & |
---|
| 88 | atime,form='unformatted',access='direct',status=file_stat,recl=24) |
---|
| 89 | |
---|
| 90 | |
---|
| 91 | ! Write current time to file |
---|
| 92 | !*************************** |
---|
| 93 | |
---|
| 94 | ! write(unitpartout_average) itime,numpart |
---|
| 95 | do i=1,numpart |
---|
| 96 | |
---|
| 97 | ! Take only valid particles |
---|
| 98 | !************************** |
---|
| 99 | |
---|
| 100 | if (itra1(i).eq.itime) then |
---|
| 101 | part_av_topo(i)=part_av_topo(i)/float(npart_av(i)) |
---|
| 102 | part_av_z(i)=part_av_z(i)/float(npart_av(i)) |
---|
| 103 | part_av_pv(i)=part_av_pv(i)/float(npart_av(i)) |
---|
| 104 | part_av_qv(i)=part_av_qv(i)/float(npart_av(i)) |
---|
| 105 | part_av_tt(i)=part_av_tt(i)/float(npart_av(i)) |
---|
| 106 | part_av_uu(i)=part_av_uu(i)/float(npart_av(i)) |
---|
| 107 | part_av_vv(i)=part_av_vv(i)/float(npart_av(i)) |
---|
| 108 | part_av_rho(i)=part_av_rho(i)/float(npart_av(i)) |
---|
| 109 | part_av_tro(i)=part_av_tro(i)/float(npart_av(i)) |
---|
| 110 | part_av_hmix(i)=part_av_hmix(i)/float(npart_av(i)) |
---|
| 111 | part_av_energy(i)=part_av_energy(i)/float(npart_av(i)) |
---|
| 112 | |
---|
| 113 | part_av_cartx(i)=part_av_cartx(i)/float(npart_av(i)) |
---|
| 114 | part_av_carty(i)=part_av_carty(i)/float(npart_av(i)) |
---|
| 115 | part_av_cartz(i)=part_av_cartz(i)/float(npart_av(i)) |
---|
| 116 | |
---|
| 117 | ! Project Cartesian coordinates back onto Earth's surface |
---|
| 118 | ! ####################################################### |
---|
| 119 | xlon=atan2(part_av_cartx(i),-1.*part_av_carty(i)) |
---|
| 120 | ylat=atan2(part_av_cartz(i),sqrt(part_av_cartx(i)*part_av_cartx(i)+ & |
---|
| 121 | part_av_carty(i)*part_av_carty(i))) |
---|
| 122 | xlon=xlon/pi180 |
---|
| 123 | ylat=ylat/pi180 |
---|
| 124 | |
---|
| 125 | |
---|
| 126 | ! Convert all data to integer*2 variables (from -32768 to 32767) for compressed output |
---|
| 127 | !************************************************************************************* |
---|
| 128 | |
---|
| 129 | if (xlon.gt.180.) xlon=xlon-360. |
---|
| 130 | if (xlon.lt.-180.) xlon=xlon+360. |
---|
| 131 | ishort_xlon(i)=nint(xlon*180.) |
---|
| 132 | ishort_ylat(i)=nint(ylat*360.) |
---|
| 133 | |
---|
| 134 | zlim=(part_av_z(i)*2.-32000.) |
---|
| 135 | zlim=min(zlim,32766.) |
---|
| 136 | zlim=max(zlim,-32766.) |
---|
| 137 | ishort_z(i)=nint(zlim) |
---|
| 138 | |
---|
| 139 | zlim=(part_av_topo(i)*2.-32000.) |
---|
| 140 | zlim=min(zlim,32766.) |
---|
| 141 | zlim=max(zlim,-32766.) |
---|
| 142 | ishort_topo(i)=nint(zlim) |
---|
| 143 | |
---|
| 144 | zlim=(part_av_tro(i)*2.-32000.) |
---|
| 145 | zlim=min(zlim,32766.) |
---|
| 146 | zlim=max(zlim,-32766.) |
---|
| 147 | ishort_tro(i)=nint(zlim) |
---|
| 148 | |
---|
| 149 | zlim=(part_av_hmix(i)*2.-32000.) |
---|
| 150 | zlim=min(zlim,32766.) |
---|
| 151 | zlim=max(zlim,-32766.) |
---|
| 152 | ishort_hmix(i)=nint(zlim) |
---|
| 153 | |
---|
| 154 | zlim=(part_av_rho(i)*20000.-32000.) |
---|
| 155 | zlim=min(zlim,32766.) |
---|
| 156 | zlim=max(zlim,-32766.) |
---|
| 157 | ishort_rho(i)=nint(zlim) |
---|
| 158 | |
---|
| 159 | zlim=(part_av_qv(i)*1000000.-32000.) |
---|
| 160 | zlim=min(zlim,32766.) |
---|
| 161 | zlim=max(zlim,-32766.) |
---|
| 162 | ishort_qv(i)=nint(zlim) |
---|
| 163 | |
---|
| 164 | zlim=(part_av_pv(i)*100.) |
---|
| 165 | zlim=min(zlim,32766.) |
---|
| 166 | zlim=max(zlim,-32766.) |
---|
| 167 | ishort_pv(i)=nint(zlim) |
---|
| 168 | |
---|
| 169 | zlim=((part_av_tt(i)-273.15))*300. |
---|
| 170 | zlim=min(zlim,32766.) |
---|
| 171 | zlim=max(zlim,-32766.) |
---|
| 172 | ishort_tt(i)=nint(zlim) |
---|
| 173 | |
---|
| 174 | zlim=(part_av_uu(i)*200.) |
---|
| 175 | zlim=min(zlim,32766.) |
---|
| 176 | zlim=max(zlim,-32766.) |
---|
| 177 | ishort_uu(i)=nint(zlim) |
---|
| 178 | |
---|
| 179 | zlim=(part_av_vv(i)*200.) |
---|
| 180 | zlim=min(zlim,32766.) |
---|
| 181 | zlim=max(zlim,-32766.) |
---|
| 182 | ishort_vv(i)=nint(zlim) |
---|
| 183 | |
---|
| 184 | zlim=(part_av_energy(i)-300000.)/30. |
---|
| 185 | zlim=min(zlim,32766.) |
---|
| 186 | zlim=max(zlim,-32766.) |
---|
| 187 | ishort_energy(i)=nint(zlim) |
---|
| 188 | |
---|
| 189 | |
---|
| 190 | ! Turn on for readable test output |
---|
| 191 | !********************************* |
---|
| 192 | |
---|
| 193 | ! write(119,*) itime,i,xlon,ylat,part_av_z(i),part_av_topo(i),part_av_tro(i), & |
---|
| 194 | ! part_av_hmix(i),part_av_rho(i),part_av_qv(i),part_av_pv(i),part_av_tt(i), & |
---|
| 195 | ! ishort_uu(i),ishort_vv(i) |
---|
| 196 | endif |
---|
| 197 | |
---|
| 198 | ! Re-initialize averages, and set number of steps to zero |
---|
| 199 | npart_av(i)=0 |
---|
| 200 | part_av_topo(i)=0. |
---|
| 201 | part_av_z(i)=0. |
---|
| 202 | part_av_cartx(i)=0. |
---|
| 203 | part_av_carty(i)=0. |
---|
| 204 | part_av_cartz(i)=0. |
---|
| 205 | part_av_pv(i)=0. |
---|
| 206 | part_av_qv(i)=0. |
---|
| 207 | part_av_tt(i)=0. |
---|
| 208 | part_av_uu(i)=0. |
---|
| 209 | part_av_vv(i)=0. |
---|
| 210 | part_av_rho(i)=0. |
---|
| 211 | part_av_tro(i)=0. |
---|
| 212 | part_av_hmix(i)=0. |
---|
| 213 | part_av_energy(i)=0. |
---|
| 214 | |
---|
| 215 | end do |
---|
| 216 | |
---|
| 217 | ! Write the output |
---|
| 218 | !***************** |
---|
| 219 | |
---|
| 220 | do i=1,numpart |
---|
| 221 | if (itra1(i).eq.itime) then |
---|
| 222 | write(unitpartout_average,rec=irec+i) ishort_xlon(i),ishort_ylat(i),ishort_z(i), & |
---|
| 223 | ishort_topo(i),ishort_tro(i),ishort_hmix(i),ishort_rho(i),ishort_qv(i),ishort_pv(i), & |
---|
| 224 | ishort_tt(i),ishort_uu(i),ishort_vv(i) ! ,ishort_energy(i) |
---|
| 225 | endif |
---|
| 226 | enddo |
---|
| 227 | |
---|
| 228 | close(unitpartout_average) |
---|
| 229 | |
---|
| 230 | end subroutine partoutput_average |
---|