source: flexpart.git/src/ohreaction.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 5 years ago

move license from headers to a different file

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