source: flexpart.git/src/ohreaction.f90 @ 02095e3

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 02095e3 was 78e62dc, checked in by flexpart <>, 9 years ago

New OH parameter in SPECIES files (now 3 instead of 2). New path to OH binariy files.

  • Property mode set to 100644
File size: 6.9 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 ohreaction(itime,ltsample,loutnext)
23  !                     i      i        i
24  !*****************************************************************************
25  !                                                                            *
26  !                                                                            *
27  !    Author: R.L. Thompson                                                   *
28  !                                                                            *
29  !    Nov 2014                                                                *
30  !                                                                            *
31  !                                                                            *
32  !*****************************************************************************
33  ! Variables:                                                                 *
34  ! ix,jy                indices of output grid cell for each particle         *
35  ! itime [s]            actual simulation time [s]                            *
36  ! jpart                particle index                                        *
37  ! ldeltat [s]          interval since radioactive decay was computed         *
38  ! loutnext [s]         time for which gridded deposition is next output      *
39  ! loutstep [s]         interval at which gridded deposition is output        *
40  ! oh_average [molecule/cm^3] OH Concentration                                *
41  ! ltsample [s]         interval over which mass is deposited                 *
42  !                                                                            *
43  !*****************************************************************************
44
45  use oh_mod
46  use par_mod
47  use com_mod
48
49  implicit none
50
51  integer :: jpart,itime,ltsample,loutnext,ldeltat,j,k,ix,jy!,ijx,jjy
52  integer :: ngrid,interp_time,n,m,h,indz,i!,ia,il
53  integer :: jjjjmmdd,hhmmss,OHx,OHy,OHz
54  real, dimension(nzOH) :: altOHtop
55  real :: xlon,ylat
56  real :: xtn,ytn
57  real :: restmass,ohreacted,oh_average
58  real :: ohrate,temp
59  real, parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
60  real(kind=dp) :: jul
61
62  ! Compute interval since radioactive decay of deposited mass was computed
63  !************************************************************************
64
65  if (itime.le.loutnext) then
66    ldeltat=itime-(loutnext-loutstep)
67  else                                  ! first half of next interval
68    ldeltat=itime-loutnext
69  endif
70
71  jul=bdate+real(itime,kind=dp)/86400.
72  call caldate(jul,jjjjmmdd,hhmmss)
73  m=(jjjjmmdd-(jjjjmmdd/10000)*10000)/100
74  h=hhmmss/10000
75
76  ! Loop over particles
77  !*****************************************
78
79  do jpart=1,numpart
80
81    ! Determine which nesting level to be used
82    ngrid=0
83    do j=numbnests,1,-1
84      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
85           (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
86        ngrid=j
87        goto 23
88      endif
89    end do
9023  continue
91
92    ! Determine nested grid coordinates
93    if (ngrid.gt.0) then
94      xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid)
95      ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid)
96      ix=int(xtn)
97      jy=int(ytn)
98    else
99      ix=int(xtra1(jpart))
100      jy=int(ytra1(jpart))
101    endif
102
103    interp_time=nint(itime-0.5*ltsample)
104    n=2
105    if(abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) n=1
106
107    do i=2,nz
108      if (height(i).gt.ztra1(jpart)) then
109        indz=i-1
110        goto 6
111      endif
112    end do
1136   continue
114
115    ! Get OH from nearest grid-cell and specific month
116    !*************************************************
117
118    ! world coordinates
119    xlon=xtra1(jpart)*dx+xlon0
120    if (xlon.gt.180) then
121       xlon=xlon-360
122    endif
123    ylat=ytra1(jpart)*dy+ylat0
124
125    ! get position in the OH field
126    OHx=minloc(abs(lonOH-xlon),dim=1,mask=abs(lonOH-xlon).eq.minval(abs(lonOH-xlon)))
127    OHy=minloc(abs(latOH-ylat),dim=1,mask=abs(latOH-ylat).eq.minval(abs(latOH-ylat)))
128
129    ! get the level of the OH field for the particle
130    ! ztra1 is the z-coord of the trajectory above model orography in metres
131    ! altOH is the height of the centre of the level in the OH field above orography
132    do i=2,nzOH
133      altOHtop(i-1)=altOH(i)+0.5*(altOH(i)-altOH(i-1))
134    end do
135    altOHtop(nzOH)=altOH(nzOH)+0.5*(altOH(nzOH)-altOH(nzOH-1))
136    OHz=minloc(abs(altOHtop-ztra1(jpart)),dim=1,mask=abs(altOHtop-ztra1(jpart))&
137            &.eq.minval(abs(altOHtop-ztra1(jpart))))
138
139    ! Interpolate between hourly OH fields to current time
140    !*****************************************************
141
142    oh_average=OH_hourly(OHx,OHy,OHz,1)+&
143               &(OH_hourly(OHx,OHy,OHz,2)-OH_hourly(OHx,OHy,OHz,1))*&
144               &(itime-memOHtime(1))/(memOHtime(2)-memOHtime(1))
145
146    if (oh_average.gt.smallnum) then
147
148      ! Computation of the OH reaction
149      !**********************************************************
150
151      temp=tt(ix,jy,indz,n)
152
153      do k=1,nspec                               
154        if (ohcconst(k).gt.0.) then
155          ohrate=ohcconst(k)*temp**ohnconst(k)*exp(-ohdconst(k)/temp)*oh_average
156          ! new particle mass
157          restmass = xmass1(jpart,k)*exp(-1*ohrate*abs(ltsample))
158          if (restmass .gt. smallnum) then
159            xmass1(jpart,k)=restmass
160          else
161            xmass1(jpart,k)=0.
162          endif
163          ohreacted=xmass1(jpart,k)*(1-exp(-1*ohrate*abs(ltsample)))
164        else
165          ohreacted=0.
166        endif
167      end do
168
169    endif  ! oh_average.gt.smallnum
170
171  end do  !continue loop over all particles
172
173
174end subroutine ohreaction
175
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG