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 readOHfield |
---|
23 | |
---|
24 | !***************************************************************************** |
---|
25 | ! * |
---|
26 | ! Reads the OH field into memory * |
---|
27 | ! * |
---|
28 | ! AUTHOR: R.L. Thompson, Nov 2014 * |
---|
29 | ! * |
---|
30 | !***************************************************************************** |
---|
31 | ! * |
---|
32 | ! Variables: * |
---|
33 | ! * |
---|
34 | ! path(numpath) contains the path names * |
---|
35 | ! lonOH(nxOH) longitude of OH fields * |
---|
36 | ! latOH(nyOH) latitude of OH fields * |
---|
37 | ! altOH(nzOH) altitude of OH fields * |
---|
38 | ! etaOH(nzOH) eta-levels of OH fields * |
---|
39 | ! OH_field(nxOH,nyOH,nzOH,m) OH concentration (molecules/cm3) * |
---|
40 | ! * |
---|
41 | ! * |
---|
42 | !***************************************************************************** |
---|
43 | |
---|
44 | use oh_mod |
---|
45 | use par_mod |
---|
46 | use com_mod |
---|
47 | |
---|
48 | implicit none |
---|
49 | |
---|
50 | include 'netcdf.inc' |
---|
51 | |
---|
52 | character(len=150) :: thefile |
---|
53 | character(len=2) :: mm |
---|
54 | integer :: nid,ierr,xid,yid,zid,vid,m |
---|
55 | real, dimension(:), allocatable :: etaOH |
---|
56 | |
---|
57 | ! real, parameter :: gasct=8.314 ! gas constant |
---|
58 | ! real, parameter :: mct=0.02894 ! kg mol-1 |
---|
59 | ! real, parameter :: g=9.80665 ! m s-2 |
---|
60 | ! real, parameter :: lrate=0.0065 ! K m-1 |
---|
61 | real, parameter :: scalehgt=7000. ! scale height in metres |
---|
62 | |
---|
63 | ! Read OH fields and level heights |
---|
64 | !******************************** |
---|
65 | |
---|
66 | do m=1,12 |
---|
67 | |
---|
68 | ! open netcdf file |
---|
69 | write(mm,fmt='(i2.2)') m |
---|
70 | ! thefile=trim(path(1))//'OH_FIELDS/'//'geos-chem.OH.2005'//mm//'01.nc' |
---|
71 | thefile=trim(ohfields_path)//'OH_FIELDS/'//'geos-chem.OH.2005'//mm//'01.nc' |
---|
72 | ierr=nf_open(trim(thefile),NF_NOWRITE,nid) |
---|
73 | if(ierr.ne.0) then |
---|
74 | write (*,*) 'The OH field at: '//thefile// ' could not be opened' |
---|
75 | write (*,*) 'please copy the OH fields there, or change the path in the' |
---|
76 | write (*,*) 'COMMAND namelist!' |
---|
77 | write(*,*) nf_strerror(ierr) |
---|
78 | stop |
---|
79 | endif |
---|
80 | |
---|
81 | ! inquire about variables |
---|
82 | ierr=nf_inq_dimid(nid,'Lon-000',xid) |
---|
83 | if(ierr.ne.0) then |
---|
84 | write(*,*) nf_strerror(ierr) |
---|
85 | stop |
---|
86 | endif |
---|
87 | ierr=nf_inq_dimid(nid,'Lat-000',yid) |
---|
88 | if(ierr.ne.0) then |
---|
89 | write(*,*) nf_strerror(ierr) |
---|
90 | stop |
---|
91 | endif |
---|
92 | ierr=nf_inq_dimid(nid,'Alt-000',zid) |
---|
93 | if(ierr.ne.0) then |
---|
94 | write(*,*) nf_strerror(ierr) |
---|
95 | stop |
---|
96 | endif |
---|
97 | |
---|
98 | if(m.eq.1) then |
---|
99 | |
---|
100 | ! read dimension sizes |
---|
101 | ierr=nf_inq_dimlen(nid,xid,nxOH) |
---|
102 | if(ierr.ne.0) then |
---|
103 | write(*,*) nf_strerror(ierr) |
---|
104 | stop |
---|
105 | endif |
---|
106 | ierr=nf_inq_dimlen(nid,yid,nyOH) |
---|
107 | if(ierr.ne.0) then |
---|
108 | write(*,*) nf_strerror(ierr) |
---|
109 | stop |
---|
110 | endif |
---|
111 | ierr=nf_inq_dimlen(nid,zid,nzOH) |
---|
112 | if(ierr.ne.0) then |
---|
113 | write(*,*) nf_strerror(ierr) |
---|
114 | stop |
---|
115 | endif |
---|
116 | |
---|
117 | ! allocate variables |
---|
118 | allocate(lonOH(nxOH)) |
---|
119 | allocate(latOH(nyOH)) |
---|
120 | allocate(etaOH(nzOH)) |
---|
121 | allocate(altOH(nzOH)) |
---|
122 | allocate(OH_field(nxOH,nyOH,nzOH,12)) |
---|
123 | allocate(OH_hourly(nxOH,nyOH,nzOH,2)) |
---|
124 | |
---|
125 | ! read dimension variables |
---|
126 | ierr=nf_inq_varid(nid,'LON',xid) |
---|
127 | ierr=nf_get_var_real(nid,xid,lonOH) |
---|
128 | if(ierr.ne.0) then |
---|
129 | write(*,*) nf_strerror(ierr) |
---|
130 | stop |
---|
131 | endif |
---|
132 | ierr=nf_inq_varid(nid,'LAT',yid) |
---|
133 | ierr=nf_get_var_real(nid,yid,latOH) |
---|
134 | if(ierr.ne.0) then |
---|
135 | write(*,*) nf_strerror(ierr) |
---|
136 | stop |
---|
137 | endif |
---|
138 | ierr=nf_inq_varid(nid,'ETAC',zid) |
---|
139 | ierr=nf_get_var_real(nid,zid,etaOH) |
---|
140 | if(ierr.ne.0) then |
---|
141 | write(*,*) nf_strerror(ierr) |
---|
142 | stop |
---|
143 | endif |
---|
144 | |
---|
145 | ! convert eta-level to altitude (assume surface pressure of 1010 hPa) |
---|
146 | altOH=log(1010./(etaOH*1010.))*scalehgt |
---|
147 | |
---|
148 | endif ! m.eq.1 |
---|
149 | |
---|
150 | ! read OH_field |
---|
151 | ierr=nf_inq_varid(nid,'CHEM-L_S__OH',vid) |
---|
152 | ierr=nf_get_var_real(nid,vid,OH_field(:,:,:,m)) |
---|
153 | if(ierr.ne.0) then |
---|
154 | write(*,*) nf_strerror(ierr) |
---|
155 | stop |
---|
156 | endif |
---|
157 | |
---|
158 | ierr=nf_close(nid) |
---|
159 | |
---|
160 | end do |
---|
161 | |
---|
162 | deallocate(etaOH) |
---|
163 | |
---|
164 | ! Read J(O1D) photolysis rates |
---|
165 | !******************************** |
---|
166 | |
---|
167 | ! open netcdf file |
---|
168 | ! thefile=trim(path(1))//'OH_FIELDS/jrate_average.nc' |
---|
169 | thefile=trim(ohfields_path)//'OH_FIELDS/jrate_average.nc' |
---|
170 | ierr=nf_open(trim(thefile),NF_NOWRITE,nid) |
---|
171 | if(ierr.ne.0) then |
---|
172 | write(*,*) nf_strerror(ierr) |
---|
173 | stop |
---|
174 | endif |
---|
175 | |
---|
176 | ! read dimension variables |
---|
177 | ierr=nf_inq_varid(nid,'longitude',xid) |
---|
178 | ierr=nf_get_var_real(nid,xid,lonjr) |
---|
179 | if(ierr.ne.0) then |
---|
180 | write(*,*) nf_strerror(ierr) |
---|
181 | stop |
---|
182 | endif |
---|
183 | ierr=nf_inq_varid(nid,'latitude',yid) |
---|
184 | ierr=nf_get_var_real(nid,yid,latjr) |
---|
185 | if(ierr.ne.0) then |
---|
186 | write(*,*) nf_strerror(ierr) |
---|
187 | stop |
---|
188 | endif |
---|
189 | |
---|
190 | ! read jrate_average |
---|
191 | ierr=nf_inq_varid(nid,'jrate',vid) |
---|
192 | ierr=nf_get_var_real(nid,vid,jrate_average) |
---|
193 | if(ierr.ne.0) then |
---|
194 | write(*,*) nf_strerror(ierr) |
---|
195 | stop |
---|
196 | endif |
---|
197 | |
---|
198 | ierr=nf_close(nid) |
---|
199 | |
---|
200 | return |
---|
201 | |
---|
202 | end subroutine readOHfield |
---|
203 | |
---|