source: flexpart.git/src/readOHfield.f90 @ f9ce123

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since f9ce123 was 8a65cb0, checked in by Espen Sollum ATMOS <espen@…>, 9 years ago

Added code, makefile for dev branch

  • Property mode set to 100644
File size: 6.5 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 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
49  implicit none
50
51  include 'netcdf.inc'
52
53  character(len=150) :: thefile
54  character(len=2) :: mm
55  integer :: nid,ierr,xid,yid,zid,vid,m
56  real, dimension(:), allocatable :: etaOH
57
58!  real, parameter :: gasct=8.314   ! gas constant
59!  real, parameter :: mct=0.02894   ! kg mol-1
60!  real, parameter :: g=9.80665     ! m s-2
61!  real, parameter :: lrate=0.0065  ! K m-1
62  real, parameter :: scalehgt=7000. ! scale height in metres
63
64  ! Read OH fields and level heights
65  !********************************
66
67  do m=1,12
68 
69    ! open netcdf file
70    write(mm,fmt='(i2.2)') m
71    thefile=trim(path(1))//'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(*,*) nf_strerror(ierr)
75      stop
76    endif
77
78    ! inquire about variables
79    ierr=nf_inq_dimid(nid,'Lon-000',xid)
80    if(ierr.ne.0) then
81      write(*,*) nf_strerror(ierr)
82      stop
83    endif
84    ierr=nf_inq_dimid(nid,'Lat-000',yid)
85    if(ierr.ne.0) then
86      write(*,*) nf_strerror(ierr)
87      stop
88    endif
89    ierr=nf_inq_dimid(nid,'Alt-000',zid)
90    if(ierr.ne.0) then
91      write(*,*) nf_strerror(ierr)
92      stop
93    endif
94
95    if(m.eq.1) then
96
97      ! read dimension sizes
98      ierr=nf_inq_dimlen(nid,xid,nxOH)
99      if(ierr.ne.0) then
100        write(*,*) nf_strerror(ierr)
101        stop
102      endif
103      ierr=nf_inq_dimlen(nid,yid,nyOH)
104      if(ierr.ne.0) then
105        write(*,*) nf_strerror(ierr)
106        stop
107      endif
108      ierr=nf_inq_dimlen(nid,zid,nzOH)
109      if(ierr.ne.0) then
110        write(*,*) nf_strerror(ierr)
111        stop
112      endif 
113
114      ! allocate variables
115      allocate(lonOH(nxOH))
116      allocate(latOH(nyOH))
117      allocate(etaOH(nzOH))
118      allocate(altOH(nzOH))
119      allocate(OH_field(nxOH,nyOH,nzOH,12))
120      allocate(OH_hourly(nxOH,nyOH,nzOH,2))
121
122      ! read dimension variables
123      ierr=nf_inq_varid(nid,'LON',xid)
124      ierr=nf_get_var_real(nid,xid,lonOH)
125      if(ierr.ne.0) then
126        write(*,*) nf_strerror(ierr)
127        stop
128      endif
129      ierr=nf_inq_varid(nid,'LAT',yid)
130      ierr=nf_get_var_real(nid,yid,latOH)
131      if(ierr.ne.0) then
132        write(*,*) nf_strerror(ierr)
133        stop
134      endif
135      ierr=nf_inq_varid(nid,'ETAC',zid)
136      ierr=nf_get_var_real(nid,zid,etaOH)
137      if(ierr.ne.0) then
138        write(*,*) nf_strerror(ierr)
139        stop
140      endif
141
142      ! convert eta-level to altitude (assume surface pressure of 1010 hPa)
143      altOH=log(1010./(etaOH*1010.))*scalehgt
144
145    endif ! m.eq.1
146
147    ! read OH_field
148    ierr=nf_inq_varid(nid,'CHEM-L_S__OH',vid)
149    ierr=nf_get_var_real(nid,vid,OH_field(:,:,:,m))
150    if(ierr.ne.0) then
151      write(*,*) nf_strerror(ierr)
152      stop
153    endif
154
155    ierr=nf_close(nid)
156
157  end do
158 
159  deallocate(etaOH)
160
161  ! Read J(O1D) photolysis rates
162  !******************************** 
163
164  ! open netcdf file
165  thefile=trim(path(1))//'OH_FIELDS/jrate_average.nc'
166  ierr=nf_open(trim(thefile),NF_NOWRITE,nid)
167  if(ierr.ne.0) then
168    write(*,*) nf_strerror(ierr)
169    stop
170  endif
171
172  ! read dimension variables
173  ierr=nf_inq_varid(nid,'longitude',xid)
174  ierr=nf_get_var_real(nid,xid,lonjr)
175  if(ierr.ne.0) then
176    write(*,*) nf_strerror(ierr)
177    stop
178  endif
179  ierr=nf_inq_varid(nid,'latitude',yid)
180  ierr=nf_get_var_real(nid,yid,latjr)
181  if(ierr.ne.0) then
182    write(*,*) nf_strerror(ierr)
183    stop
184  endif
185
186  ! read jrate_average
187  ierr=nf_inq_varid(nid,'jrate',vid)
188  ierr=nf_get_var_real(nid,vid,jrate_average)
189  if(ierr.ne.0) then
190    write(*,*) nf_strerror(ierr)
191    stop
192  endif
193
194  ierr=nf_close(nid)
195
196  return
197
198end subroutine readOHfield
199
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG