source: flexpart.git/src/readOHfield.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@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 3.0 KB
Line 
1subroutine readOHfield
2
3!*****************************************************************************
4!                                                                            *
5! Reads the OH field into memory                                             *
6!                                                                            *
7! AUTHOR: R.L. Thompson, Nov 2014                                            *
8!                                                                            *
9! UPDATES:                                                                   *
10!   03/2018 SEC: Converted original netCDF files to binary format            *
11!*****************************************************************************
12!                                                                            *
13! Variables:                                                                 *
14!                                                                            *
15! path(numpath)              contains the path names                         *
16! lonOH(nxOH)                longitude of OH fields                          *
17! latOH(nyOH)                latitude of OH fields                           *
18! altOH(nzOH)                altitude of OH fields                           *
19! etaOH(nzOH)                eta-levels of OH fields                         *
20! OH_field(nxOH,nyOH,nzOH,m) OH concentration (molecules/cm3)                *
21!                                                                            *
22!                                                                            *
23!*****************************************************************************
24
25  use oh_mod
26  use par_mod
27  use com_mod
28
29  implicit none
30
31  integer :: i,j,k,l,ierr
32  real, dimension(:), allocatable :: etaOH
33
34!  real, parameter :: gasct=8.314   ! gas constant
35!  real, parameter :: mct=0.02894   ! kg mol-1
36!  real, parameter :: g=9.80665     ! m s-2
37!  real, parameter :: lrate=0.0065  ! K m-1
38  real, parameter :: scalehgt=7000. ! scale height in metres
39
40
41  open(unitOH,file=trim(ohfields_path) &
42       //'OH_FIELDS/OH_variables.bin',status='old', &
43       form='UNFORMATTED', iostat=ierr, convert='little_endian')
44
45  if(ierr.ne.0) then
46    write(*,*) 'Cannot read binary OH fields in ',trim(ohfields_path)//'OH_FIELDS/OH_variables.bin'
47    stop
48  endif
49
50  read(unitOH) nxOH
51  read(unitOH) nyOH
52  read(unitOH) nzOH
53  write(*,*) nxOH,nyOH,nzOH
54
55! allocate variables
56  allocate(lonOH(nxOH))
57  allocate(latOH(nyOH))
58  allocate(etaOH(nzOH))
59  allocate(altOH(nzOH))
60  allocate(OH_field(nxOH,nyOH,nzOH,12))
61  allocate(OH_hourly(nxOH,nyOH,nzOH,2))
62
63  read(unitOH) (lonjr(i),i=1,360)
64  read(unitOH) (latjr(i),i=1,180)
65  read(unitOH) (((jrate_average(i,j,k),i=1,360),j=1,180),k=1,12)
66  read(unitOH) (lonOH(i),i=1,nxOH)
67  read(unitOH) (latOH(i),i=1,nyOH)
68  read(unitOH) (lonOH(i),i=1,nxOH)
69
70  read(unitOH) (altOH(i),i=1,nzOH)
71  read(unitOH) ((((OH_field(i,j,k,l),i=1,nxOH),j=1,nyOH),k=1,nzOH),l=1,12)
72  read(unitOH) ((((OH_hourly(i,j,k,l),i=1,nxOH),j=1,nyOH),k=1,nzOH),l=1,2)
73
74end subroutine readOHfield
75
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG