source: flexpart.git/src/calcfluxes.f90 @ a7167e4

dev
Last change on this file since a7167e4 was a7167e4, checked in by Sabine <sabine.eckhardt@…>, 4 years ago

bugfix: netcdffluxfield init, flux through upper, southern and western boundary

  • Property mode set to 100644
File size: 6.3 KB
Line 
1subroutine calcfluxes(nage,jpart,xold,yold,zold)
2  !                       i     i    i    i    i
3  !*****************************************************************************
4  !                                                                            *
5  !     Calculation of the gross fluxes across horizontal, eastward and        *
6  !     northward facing surfaces. The routine calculates the mass flux        *
7  !     due to the motion of only one particle. The fluxes of subsequent calls *
8  !     to this subroutine are accumulated until the next output is due.       *
9  !     Upon output, flux fields are re-set to zero in subroutine fluxoutput.f.*
10  !                                                                            *
11  !     Author: A. Stohl                                                       *
12  !                                                                            *
13  !     04 April 2000                                                          *
14  !                                                                            *
15  !*****************************************************************************
16  !                                                                            *
17  ! Variables:                                                                 *
18  !                                                                            *
19  ! nage                  Age class of the particle considered                 *
20  ! jpart                 Index of the particle considered                     *
21  ! xold,yold,zold        "Memorized" old positions of the particle            *
22  !                                                                            *
23  !*****************************************************************************
24
25  use flux_mod
26  use outg_mod
27  use par_mod
28  use com_mod
29
30  implicit none
31
32  integer :: jpart,nage,ixave,jyave,kz,kzave,kp
33  integer :: k,k1,k2,ix,ix1,ix2,ixs,jy,jy1,jy2
34  real :: xold,yold,zold,xmean,ymean
35
36
37  ! Determine average positions
38  !****************************
39
40  if ((ioutputforeachrelease.eq.1).and.(mdomainfill.eq.0)) then
41     kp=npoint(jpart)
42  else
43     kp=1
44  endif
45
46  xmean=(xold+xtra1(jpart))/2.
47  ymean=(yold+ytra1(jpart))/2.
48
49  ixave=int((xmean*dx+xoutshift)/dxout)
50  jyave=int((ymean*dy+youtshift)/dyout)
51  do kz=1,numzgrid                ! determine height of cell
52    if (outheight(kz).gt.ztra1(jpart)) goto 16
53  end do
5416   kzave=kz
55
56
57  ! Determine vertical fluxes
58  !**************************
59
60  if ((ixave.ge.0).and.(jyave.ge.0).and.(ixave.le.numxgrid-1).and. &
61       (jyave.le.numygrid-1)) then
62    do kz=1,numzgrid                ! determine height of cell
63!      if (outheighthalf(kz).gt.zold) goto 11
64      if (outheight(kz).gt.zold) goto 11 !sec, use upper layer instead of mid layer
65    end do
6611   k1=min(numzgrid,kz)
67    do kz=1,numzgrid                ! determine height of cell
68!      if (outheighthalf(kz).gt.ztra1(jpart)) goto 21
69      if (outheight(kz).gt.ztra1(jpart)) goto 21
70    end do
7121   k2=min(numzgrid,kz)
72
73    do k=1,nspec
74      do kz=k1,k2-1
75        flux(5,ixave,jyave,kz,k,kp,nage)= &
76             flux(5,ixave,jyave,kz,k,kp,nage)+ &
77             xmass1(jpart,k)
78      end do
79      do kz=k2,k1-1
80        flux(6,ixave,jyave,kz,k,kp,nage)= &
81             flux(6,ixave,jyave,kz,k,kp,nage)+ &
82             xmass1(jpart,k)
83      end do
84    end do
85  endif
86
87
88  ! Determine west-east fluxes (fluxw) and east-west fluxes (fluxe)
89  !****************************************************************
90
91  if ((kzave.le.numzgrid).and.(jyave.ge.0).and. &
92       (jyave.le.numygrid-1)) then
93
94  ! 1) Particle does not cross domain boundary
95
96    if (abs(xold-xtra1(jpart)).lt.real(nx)/2.) then
97      ix1=int((xold*dx+xoutshift)/dxout)              ! flux throught the western boundary (sec)
98      ix2=int((xtra1(jpart)*dx+xoutshift)/dxout)      ! flux throught the western boundary (sec)
99!      ix1=int((xold*dx+xoutshift)/dxout+0.5)
100!      ix2=int((xtra1(jpart)*dx+xoutshift)/dxout+0.5)
101      do k=1,nspec
102        do ix=ix1,ix2-1
103          if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
104            flux(1,ix,jyave,kzave,k,kp,nage)= &
105                 flux(1,ix,jyave,kzave,k,kp,nage) &
106                 +xmass1(jpart,k)
107          endif
108        end do
109        do ix=ix2,ix1-1
110          if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
111            flux(2,ix,jyave,kzave,k,kp,nage)= &
112                 flux(2,ix,jyave,kzave,k,kp,nage) &
113                 +xmass1(jpart,k)
114          endif
115        end do
116      end do
117
118  ! 2) Particle crosses domain boundary: use cyclic boundary condition
119  !    and attribute flux to easternmost grid row only (approximation valid
120  !    for relatively slow motions compared to output grid cell size)
121
122    else
123      ixs=int(((real(nxmin1)-1.e5)*dx+xoutshift)/dxout)
124      if ((ixs.ge.0).and.(ixs.le.numxgrid-1)) then
125        if (xold.gt.xtra1(jpart)) then       ! west-east flux
126          do k=1,nspec
127            flux(1,ixs,jyave,kzave,k,kp,nage)= &
128                 flux(1,ixs,jyave,kzave,k,kp,nage) &
129                 +xmass1(jpart,k)
130          end do
131        else                                 ! east-west flux
132          do k=1,nspec
133            flux(2,ixs,jyave,kzave,k,kp,nage)= &
134                 flux(2,ixs,jyave,kzave,k,kp,nage) &
135                 +xmass1(jpart,k)
136          end do
137        endif
138      endif
139    endif
140  endif
141
142
143  ! Determine south-north fluxes (fluxs) and north-south fluxes (fluxn)
144  !********************************************************************
145
146  if ((kzave.le.numzgrid).and.(ixave.ge.0).and. &
147       (ixave.le.numxgrid-1)) then
148    jy1=int((yold*dy+youtshift)/dyout)         ! flux throught the southern boundary (sec)
149    jy2=int((ytra1(jpart)*dy+youtshift)/dyout) ! flux throught the southern boundary (sec)
150!    jy1=int((yold*dy+youtshift)/dyout+0.5)
151!    jy2=int((ytra1(jpart)*dy+youtshift)/dyout+0.5)
152
153    do k=1,nspec
154      do jy=jy1,jy2-1
155        if ((jy.ge.0).and.(jy.le.numygrid-1)) then
156          flux(3,ixave,jy,kzave,k,kp,nage)= &
157               flux(3,ixave,jy,kzave,k,kp,nage) &
158               +xmass1(jpart,k)
159        endif
160      end do
161      do jy=jy2,jy1-1
162        if ((jy.ge.0).and.(jy.le.numygrid-1)) then
163          flux(4,ixave,jy,kzave,k,kp,nage)= &
164               flux(4,ixave,jy,kzave,k,kp,nage) &
165               +xmass1(jpart,k)
166        endif
167      end do
168    end do
169  endif
170
171end subroutine calcfluxes
172
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG