[0c8c7f2] | 1 | |
---|
| 2 | subroutine partpos_average(itime,j) |
---|
| 3 | |
---|
| 4 | |
---|
| 5 | !********************************************************************** |
---|
| 6 | ! This subroutine averages particle quantities, to be used for particle |
---|
| 7 | ! dump (in partoutput.f90). Averaging is done over output interval. |
---|
| 8 | !********************************************************************** |
---|
| 9 | |
---|
| 10 | use par_mod |
---|
| 11 | use com_mod |
---|
| 12 | |
---|
| 13 | implicit none |
---|
| 14 | |
---|
| 15 | integer :: itime,j,ix,jy,ixp,jyp,indexh,m,il,ind,indz,indzp |
---|
| 16 | real :: xlon,ylat,x,y,z |
---|
| 17 | real :: dt1,dt2,dtt,ddx,ddy,rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz |
---|
| 18 | real :: topo,hm(2),hmixi,pv1(2),pvprof(2),pvi,qv1(2),qvprof(2),qvi |
---|
| 19 | real :: tt1(2),ttprof(2),tti,rho1(2),rhoprof(2),rhoi |
---|
| 20 | real :: uu1(2),uuprof(2),uui,vv1(2),vvprof(2),vvi |
---|
| 21 | real :: tr(2),tri,energy |
---|
| 22 | |
---|
| 23 | |
---|
| 24 | |
---|
| 25 | ! Some variables needed for temporal interpolation |
---|
| 26 | !************************************************* |
---|
| 27 | |
---|
| 28 | dt1=real(itime-memtime(1)) |
---|
| 29 | dt2=real(memtime(2)-itime) |
---|
| 30 | dtt=1./(dt1+dt2) |
---|
| 31 | |
---|
| 32 | xlon=xlon0+xtra1(j)*dx |
---|
| 33 | ylat=ylat0+ytra1(j)*dy |
---|
| 34 | |
---|
| 35 | !***************************************************************************** |
---|
| 36 | ! Interpolate several variables (PV, specific humidity, etc.) to particle position |
---|
| 37 | !***************************************************************************** |
---|
| 38 | |
---|
| 39 | ix=xtra1(j) |
---|
| 40 | jy=ytra1(j) |
---|
| 41 | ixp=ix+1 |
---|
| 42 | jyp=jy+1 |
---|
| 43 | ddx=xtra1(j)-real(ix) |
---|
| 44 | ddy=ytra1(j)-real(jy) |
---|
| 45 | rddx=1.-ddx |
---|
| 46 | rddy=1.-ddy |
---|
| 47 | p1=rddx*rddy |
---|
| 48 | p2=ddx*rddy |
---|
| 49 | p3=rddx*ddy |
---|
| 50 | p4=ddx*ddy |
---|
| 51 | |
---|
| 52 | ! eso: Temporary fix for particle exactly at north pole |
---|
| 53 | if (jyp >= nymax) then |
---|
| 54 | ! write(*,*) 'WARNING: conccalc.f90 jyp >= nymax' |
---|
| 55 | jyp=jyp-1 |
---|
| 56 | end if |
---|
| 57 | |
---|
| 58 | ! Topography |
---|
| 59 | !*********** |
---|
| 60 | |
---|
| 61 | topo=p1*oro(ix,jy)+p2*oro(ixp,jy)+p3*oro(ix,jyp)+p4*oro(ixp,jyp) |
---|
| 62 | |
---|
| 63 | ! Potential vorticity, specific humidity, temperature, and density |
---|
| 64 | !***************************************************************** |
---|
| 65 | |
---|
| 66 | do il=2,nz |
---|
| 67 | if (height(il).gt.ztra1(j)) then |
---|
| 68 | indz=il-1 |
---|
| 69 | indzp=il |
---|
| 70 | goto 6 |
---|
| 71 | endif |
---|
| 72 | end do |
---|
| 73 | 6 continue |
---|
| 74 | |
---|
| 75 | dz1=ztra1(j)-height(indz) |
---|
| 76 | dz2=height(indzp)-ztra1(j) |
---|
| 77 | dz=1./(dz1+dz2) |
---|
| 78 | |
---|
| 79 | |
---|
| 80 | do ind=indz,indzp |
---|
| 81 | do m=1,2 |
---|
| 82 | indexh=memind(m) |
---|
| 83 | |
---|
| 84 | ! Potential vorticity |
---|
| 85 | pv1(m)=p1*pv(ix ,jy ,ind,indexh) & |
---|
| 86 | +p2*pv(ixp,jy ,ind,indexh) & |
---|
| 87 | +p3*pv(ix ,jyp,ind,indexh) & |
---|
| 88 | +p4*pv(ixp,jyp,ind,indexh) |
---|
| 89 | ! Specific humidity |
---|
| 90 | qv1(m)=p1*qv(ix ,jy ,ind,indexh) & |
---|
| 91 | +p2*qv(ixp,jy ,ind,indexh) & |
---|
| 92 | +p3*qv(ix ,jyp,ind,indexh) & |
---|
| 93 | +p4*qv(ixp,jyp,ind,indexh) |
---|
| 94 | ! Temperature |
---|
| 95 | tt1(m)=p1*tt(ix ,jy ,ind,indexh) & |
---|
| 96 | +p2*tt(ixp,jy ,ind,indexh) & |
---|
| 97 | +p3*tt(ix ,jyp,ind,indexh) & |
---|
| 98 | +p4*tt(ixp,jyp,ind,indexh) |
---|
| 99 | ! U wind |
---|
| 100 | uu1(m)=p1*uu(ix ,jy ,ind,indexh) & |
---|
| 101 | +p2*uu(ixp,jy ,ind,indexh) & |
---|
| 102 | +p3*uu(ix ,jyp,ind,indexh) & |
---|
| 103 | +p4*uu(ixp,jyp,ind,indexh) |
---|
| 104 | ! V wind |
---|
| 105 | vv1(m)=p1*vv(ix ,jy ,ind,indexh) & |
---|
| 106 | +p2*vv(ixp,jy ,ind,indexh) & |
---|
| 107 | +p3*vv(ix ,jyp,ind,indexh) & |
---|
| 108 | +p4*vv(ixp,jyp,ind,indexh) |
---|
| 109 | ! Density |
---|
| 110 | rho1(m)=p1*rho(ix ,jy ,ind,indexh) & |
---|
| 111 | +p2*rho(ixp,jy ,ind,indexh) & |
---|
| 112 | +p3*rho(ix ,jyp,ind,indexh) & |
---|
| 113 | +p4*rho(ixp,jyp,ind,indexh) |
---|
| 114 | end do |
---|
| 115 | pvprof(ind-indz+1)=(pv1(1)*dt2+pv1(2)*dt1)*dtt |
---|
| 116 | qvprof(ind-indz+1)=(qv1(1)*dt2+qv1(2)*dt1)*dtt |
---|
| 117 | ttprof(ind-indz+1)=(tt1(1)*dt2+tt1(2)*dt1)*dtt |
---|
| 118 | uuprof(ind-indz+1)=(uu1(1)*dt2+uu1(2)*dt1)*dtt |
---|
| 119 | vvprof(ind-indz+1)=(vv1(1)*dt2+vv1(2)*dt1)*dtt |
---|
| 120 | rhoprof(ind-indz+1)=(rho1(1)*dt2+rho1(2)*dt1)*dtt |
---|
| 121 | end do |
---|
| 122 | pvi=(dz1*pvprof(2)+dz2*pvprof(1))*dz |
---|
| 123 | qvi=(dz1*qvprof(2)+dz2*qvprof(1))*dz |
---|
| 124 | tti=(dz1*ttprof(2)+dz2*ttprof(1))*dz |
---|
| 125 | uui=(dz1*uuprof(2)+dz2*uuprof(1))*dz |
---|
| 126 | vvi=(dz1*vvprof(2)+dz2*vvprof(1))*dz |
---|
| 127 | rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz |
---|
| 128 | |
---|
| 129 | ! Tropopause and PBL height |
---|
| 130 | !************************** |
---|
| 131 | |
---|
| 132 | do m=1,2 |
---|
| 133 | indexh=memind(m) |
---|
| 134 | |
---|
| 135 | ! Tropopause |
---|
| 136 | tr(m)=p1*tropopause(ix ,jy ,1,indexh) & |
---|
| 137 | + p2*tropopause(ixp,jy ,1,indexh) & |
---|
| 138 | + p3*tropopause(ix ,jyp,1,indexh) & |
---|
| 139 | + p4*tropopause(ixp,jyp,1,indexh) |
---|
| 140 | |
---|
| 141 | ! PBL height |
---|
| 142 | hm(m)=p1*hmix(ix ,jy ,1,indexh) & |
---|
| 143 | + p2*hmix(ixp,jy ,1,indexh) & |
---|
| 144 | + p3*hmix(ix ,jyp,1,indexh) & |
---|
| 145 | + p4*hmix(ixp,jyp,1,indexh) |
---|
| 146 | end do |
---|
| 147 | |
---|
| 148 | hmixi=(hm(1)*dt2+hm(2)*dt1)*dtt |
---|
| 149 | tri=(tr(1)*dt2+tr(2)*dt1)*dtt |
---|
| 150 | |
---|
| 151 | |
---|
| 152 | energy=tti*cpa+(ztra1(j)+topo)*9.81+qvi*2501000.+(uui**2+vvi**2)/2. |
---|
| 153 | |
---|
| 154 | ! Add new values to sum and increase counter by one |
---|
| 155 | !************************************************** |
---|
| 156 | |
---|
| 157 | npart_av(j)=npart_av(j)+1 |
---|
| 158 | |
---|
| 159 | ! Calculate Cartesian 3D coordinates suitable for averaging |
---|
| 160 | !********************************************************** |
---|
| 161 | |
---|
| 162 | xlon=xlon*pi180 |
---|
| 163 | ylat=ylat*pi180 |
---|
| 164 | x = cos(ylat)*sin(xlon) |
---|
| 165 | y = -1.*cos(ylat)*cos(xlon) |
---|
| 166 | z = sin(ylat) |
---|
| 167 | |
---|
| 168 | part_av_cartx(j)=part_av_cartx(j)+x |
---|
| 169 | part_av_carty(j)=part_av_carty(j)+y |
---|
| 170 | part_av_cartz(j)=part_av_cartz(j)+z |
---|
| 171 | part_av_z(j)=part_av_z(j)+ztra1(j) |
---|
| 172 | part_av_topo(j)=part_av_topo(j)+topo |
---|
| 173 | part_av_pv(j)=part_av_pv(j)+pvi |
---|
| 174 | part_av_qv(j)=part_av_qv(j)+qvi |
---|
| 175 | part_av_tt(j)=part_av_tt(j)+tti |
---|
| 176 | part_av_uu(j)=part_av_uu(j)+uui |
---|
| 177 | part_av_vv(j)=part_av_vv(j)+vvi |
---|
| 178 | part_av_rho(j)=part_av_rho(j)+rhoi |
---|
| 179 | part_av_tro(j)=part_av_tro(j)+tri |
---|
| 180 | part_av_hmix(j)=part_av_hmix(j)+hmixi |
---|
| 181 | part_av_energy(j)=part_av_energy(j)+energy |
---|
| 182 | |
---|
| 183 | |
---|
| 184 | return |
---|
| 185 | end subroutine partpos_average |
---|