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