[92fab65] | 1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
| 2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
[332fbbd] | 3 | |
---|
[8a65cb0] | 4 | subroutine gethourlyOH(itime) |
---|
| 5 | ! i |
---|
| 6 | !***************************************************************************** |
---|
| 7 | ! * |
---|
| 8 | ! * |
---|
| 9 | ! Author: R.L. Thompson * |
---|
| 10 | ! * |
---|
| 11 | ! Nov 2014 * |
---|
| 12 | ! * |
---|
| 13 | ! * |
---|
| 14 | !***************************************************************************** |
---|
| 15 | ! Variables: * |
---|
| 16 | ! * |
---|
| 17 | !***************************************************************************** |
---|
| 18 | |
---|
| 19 | use oh_mod |
---|
| 20 | use par_mod |
---|
| 21 | use com_mod |
---|
| 22 | |
---|
| 23 | implicit none |
---|
| 24 | |
---|
| 25 | integer :: itime |
---|
| 26 | integer :: ix,jy,kz,m1,m2 |
---|
| 27 | integer :: ijx,jjy |
---|
| 28 | integer :: jjjjmmdd,hhmmss |
---|
| 29 | real :: sza,jrate,photo_O1D,zenithangle |
---|
| 30 | real(kind=dp) :: jul1,jul2 |
---|
| 31 | |
---|
| 32 | ! print*, 'itime: ',itime |
---|
| 33 | ! print*, 'memOHtime(1):',memOHtime(1) |
---|
| 34 | ! print*, 'memOHtime(2):',memOHtime(2) |
---|
| 35 | |
---|
| 36 | |
---|
| 37 | ! Check hourly OH field is available for the current time step |
---|
| 38 | !************************************************************** |
---|
| 39 | |
---|
| 40 | if ((ldirect*memOHtime(1).le.ldirect*itime).and. & |
---|
| 41 | (ldirect*memOHtime(2).gt.ldirect*itime)) then |
---|
| 42 | |
---|
| 43 | ! The right OH fields are already in memory -> don't do anything |
---|
| 44 | !**************************************************************** |
---|
| 45 | |
---|
| 46 | continue |
---|
| 47 | |
---|
| 48 | else if ((ldirect*memOHtime(2).le.ldirect*itime).and. & |
---|
| 49 | (memOHtime(2).ne.0.)) then |
---|
| 50 | |
---|
| 51 | ! Current time is after 2nd OH field |
---|
| 52 | !************************************ |
---|
| 53 | |
---|
| 54 | memOHtime(1)=memOHtime(2) |
---|
| 55 | memOHtime(2)=memOHtime(1)+ldirect*3600. |
---|
| 56 | OH_hourly(:,:,:,1)=OH_hourly(:,:,:,2) |
---|
| 57 | |
---|
| 58 | ! Compute new hourly value of OH |
---|
| 59 | !********************************************************** |
---|
| 60 | |
---|
| 61 | jul2=bdate+memOHtime(2)/86400._dp ! date for next hour |
---|
| 62 | call caldate(jul2,jjjjmmdd,hhmmss) |
---|
| 63 | m2=(jjjjmmdd-(jjjjmmdd/10000)*10000)/100 |
---|
| 64 | |
---|
| 65 | ! print*, 'jul2:',jul2 |
---|
| 66 | ! print*, 'm2:',m2 |
---|
| 67 | |
---|
| 68 | do kz=1,nzOH |
---|
| 69 | do jy=1,nyOH |
---|
| 70 | do ix=1,nxOH |
---|
| 71 | ijx=minloc(abs(lonjr-lonOH(ix)),dim=1,mask=abs(lonjr-lonOH(ix)).eq.minval(abs(lonjr-lonOH(ix)))) |
---|
| 72 | jjy=minloc(abs(latjr-latOH(jy)),dim=1,mask=abs(latjr-latOH(jy)).eq.minval(abs(latjr-latOH(jy)))) |
---|
| 73 | ! calculate solar zenith angle in degrees (sza) |
---|
| 74 | sza=zenithangle(latOH(jy),lonOH(ix),jul2) |
---|
| 75 | ! calculate J(O1D) (jrate) |
---|
| 76 | jrate=photo_O1D(sza) |
---|
| 77 | ! apply hourly correction to OH |
---|
| 78 | if(jrate_average(ijx,jjy,m2).gt.0.) then |
---|
| 79 | OH_hourly(ix,jy,kz,2)=OH_field(ix,jy,kz,m2)*jrate/jrate_average(ijx,jjy,m2) |
---|
| 80 | else |
---|
| 81 | OH_hourly(ix,jy,kz,2)=0. |
---|
| 82 | endif |
---|
| 83 | !! for testing !! |
---|
| 84 | ! if(jy.eq.36.and.ix.eq.36.and.kz.eq.1) then |
---|
| 85 | ! write(999,fmt='(F6.3)') jrate/jrate_average(ijx,jjy,m2) |
---|
| 86 | ! endif |
---|
| 87 | ! if(jy.eq.11.and.ix.eq.36.and.kz.eq.1) then |
---|
| 88 | ! write(998,fmt='(F6.3)') jrate/jrate_average(ijx,jjy,m2) |
---|
| 89 | ! endif |
---|
| 90 | end do |
---|
| 91 | end do |
---|
| 92 | end do |
---|
| 93 | |
---|
| 94 | else |
---|
| 95 | |
---|
| 96 | ! No OH fields in memory -> compute both hourly OH fields |
---|
| 97 | !********************************************************** |
---|
| 98 | |
---|
| 99 | jul1=bdate ! begin date of simulation (julian) |
---|
| 100 | call caldate(jul1,jjjjmmdd,hhmmss) |
---|
| 101 | m1=(jjjjmmdd-(jjjjmmdd/10000)*10000)/100 |
---|
| 102 | memOHtime(1)=0. |
---|
| 103 | |
---|
[78e62dc] | 104 | jul2=bdate+ldirect*real(1./24.,kind=dp) ! date for next hour |
---|
[8a65cb0] | 105 | call caldate(jul2,jjjjmmdd,hhmmss) |
---|
| 106 | m2=(jjjjmmdd-(jjjjmmdd/10000)*10000)/100 |
---|
| 107 | memOHtime(2)=ldirect*3600. |
---|
| 108 | |
---|
| 109 | ! print*, 'jul1:',jul1 |
---|
| 110 | ! print*, 'jul2:',jul2 |
---|
| 111 | ! print*, 'm1,m2:',m1,m2 |
---|
| 112 | |
---|
| 113 | do kz=1,nzOH |
---|
| 114 | do jy=1,nyOH |
---|
| 115 | do ix=1,nxOH |
---|
| 116 | ijx=minloc(abs(lonjr-lonOH(ix)),dim=1,mask=abs(lonjr-lonOH(ix)).eq.minval(abs(lonjr-lonOH(ix)))) |
---|
| 117 | jjy=minloc(abs(latjr-latOH(jy)),dim=1,mask=abs(latjr-latOH(jy)).eq.minval(abs(latjr-latOH(jy)))) |
---|
| 118 | ! calculate solar zenith angle in degrees (sza), beginning |
---|
| 119 | sza=zenithangle(latOH(jy),lonOH(ix),jul1) |
---|
| 120 | ! calculate J(O1D) (jrate), beginning |
---|
| 121 | jrate=photo_O1D(sza) |
---|
| 122 | ! apply hourly correction to OH |
---|
| 123 | if(jrate_average(ijx,jjy,m1).gt.0.) then |
---|
| 124 | OH_hourly(ix,jy,kz,1)=OH_field(ix,jy,kz,m1)*jrate/jrate_average(ijx,jjy,m1) |
---|
| 125 | else |
---|
| 126 | OH_hourly(ix,jy,kz,1)=0. |
---|
| 127 | endif |
---|
| 128 | ! calculate solar zenith angle in degrees (sza), after 1-hour |
---|
| 129 | sza=zenithangle(latOH(jy),lonOH(ix),jul2) |
---|
| 130 | ! calculate J(O1D) (jrate), after 1-hour |
---|
| 131 | jrate=photo_O1D(sza) |
---|
| 132 | ! apply hourly correction to OH |
---|
| 133 | if(jrate_average(ijx,jjy,m2).gt.0.) then |
---|
| 134 | OH_hourly(ix,jy,kz,2)=OH_field(ix,jy,kz,m2)*jrate/jrate_average(ijx,jjy,m2) |
---|
| 135 | else |
---|
| 136 | OH_hourly(ix,jy,kz,2)=0. |
---|
| 137 | endif |
---|
| 138 | end do |
---|
| 139 | end do |
---|
| 140 | end do |
---|
| 141 | |
---|
| 142 | endif |
---|
| 143 | |
---|
| 144 | end subroutine gethourlyOH |
---|
| 145 | |
---|