source: branches/flexpart91_hasod/src/readoutgrid.f90 @ 7

Last change on this file since 7 was 7, checked in by hasod, 11 years ago

Initial import

  • namelist input for COMMAND
  • pathnames optionally as command line argument
  • conversion utility from COMMAND to COMMAND namelist
File size: 7.1 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 readoutgrid
23
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine reads the user specifications for the output grid.        *
27  !                                                                            *
28  !     Author: A. Stohl                                                       *
29  !                                                                            *
30  !     4 June 1996                                                            *
31  !                                                                            *
32  !*****************************************************************************
33  !                                                                            *
34  ! Variables:                                                                 *
35  ! dxout,dyout          grid distance                                         *
36  ! numxgrid,numygrid,numzgrid    grid dimensions                              *
37  ! outlon0,outlat0      lower left corner of grid                             *
38  ! outheight(maxzgrid)  height levels of output grid [m]                      *
39  !                                                                            *
40  ! Constants:                                                                 *
41  ! unitoutgrid          unit connected to file OUTGRID                        *
42  !                                                                            *
43  !*****************************************************************************
44
45  use outg_mod
46  use par_mod
47  use com_mod
48
49  implicit none
50
51  integer :: i,j,stat
52  real :: outhelp,xr,xr1,yr,yr1
53  real,parameter :: eps=1.e-4
54
55
56
57  ! Open the OUTGRID file and read output grid specifications
58  !**********************************************************
59
60  open(unitoutgrid,file=path(1)(1:length(1))//'OUTGRID',status='old', &
61       err=999)
62
63
64  call skplin(5,unitoutgrid)
65
66
67  ! 1.  Read horizontal grid specifications
68  !****************************************
69
70  call skplin(3,unitoutgrid)
71  read(unitoutgrid,'(4x,f11.4)') outlon0
72  call skplin(3,unitoutgrid)
73  read(unitoutgrid,'(4x,f11.4)') outlat0
74  call skplin(3,unitoutgrid)
75  read(unitoutgrid,'(4x,i5)') numxgrid
76  call skplin(3,unitoutgrid)
77  read(unitoutgrid,'(4x,i5)') numygrid
78  call skplin(3,unitoutgrid)
79  read(unitoutgrid,'(4x,f12.5)') dxout
80  call skplin(3,unitoutgrid)
81  read(unitoutgrid,'(4x,f12.5)') dyout
82
83
84  ! Check validity of output grid (shall be within model domain)
85  !*************************************************************
86
87  xr=outlon0+real(numxgrid)*dxout
88  yr=outlat0+real(numygrid)*dyout
89  xr1=xlon0+real(nxmin1)*dx
90  yr1=ylat0+real(nymin1)*dy
91  if ((outlon0+eps.lt.xlon0).or.(outlat0+eps.lt.ylat0) &
92       .or.(xr.gt.xr1+eps).or.(yr.gt.yr1+eps)) then
93    write(*,*) outlon0,outlat0
94    write(*,*) xr1,yr1,xlon0,ylat0,xr,yr,dxout,dyout
95    write(*,*) ' #### FLEXPART MODEL ERROR! PART OF OUTPUT    ####'
96    write(*,*) ' #### GRID IS OUTSIDE MODEL DOMAIN. CHANGE    ####'
97    write(*,*) ' #### FILE OUTGRID IN DIRECTORY               ####'
98    write(*,'(a)') path(1)(1:length(1))
99    stop
100  endif
101
102  ! 2. Count Vertical levels of output grid
103  !****************************************
104  j=0
105100   j=j+1
106    do i=1,3
107      read(unitoutgrid,*,end=99)
108    end do
109    read(unitoutgrid,'(4x,f7.1)',end=99) outhelp
110    if (outhelp.eq.0.) goto 99
111    goto 100
11299   numzgrid=j-1
113
114    allocate(outheight(numzgrid) &
115         ,stat=stat)
116    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
117    allocate(outheighthalf(numzgrid) &
118         ,stat=stat)
119    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
120
121
122  rewind(unitoutgrid)
123  call skplin(29,unitoutgrid)
124
125  ! 2. Vertical levels of output grid
126  !**********************************
127
128  j=0
1291000   j=j+1
130    do i=1,3
131      read(unitoutgrid,*,end=990)
132    end do
133    read(unitoutgrid,'(4x,f7.1)',end=990) outhelp
134    if (outhelp.eq.0.) goto 99
135    outheight(j)=outhelp
136    goto 1000
137990   numzgrid=j-1
138
139
140  ! Check whether vertical levels are specified in ascending order
141  !***************************************************************
142
143  do j=2,numzgrid
144    if (outheight(j).le.outheight(j-1)) then
145    write(*,*) ' #### FLEXPART MODEL ERROR! YOUR SPECIFICATION#### '
146    write(*,*) ' #### OF OUTPUT LEVELS IS CORRUPT AT LEVEL    #### '
147    write(*,*) ' #### ',j,'                              #### '
148    write(*,*) ' #### PLEASE MAKE CHANGES IN FILE OUTGRID.    #### '
149    endif
150  end do
151
152  ! Determine the half levels, i.e. middle levels of the output grid
153  !*****************************************************************
154
155  outheighthalf(1)=outheight(1)/2.
156  do j=2,numzgrid
157    outheighthalf(j)=(outheight(j-1)+outheight(j))/2.
158  end do
159
160
161  xoutshift=xlon0-outlon0
162  youtshift=ylat0-outlat0
163  close(unitoutgrid)
164
165    allocate(oroout(0:numxgrid-1,0:numygrid-1) &
166         ,stat=stat)
167    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
168    allocate(area(0:numxgrid-1,0:numygrid-1) &
169         ,stat=stat)
170    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
171    allocate(volume(0:numxgrid-1,0:numygrid-1,numzgrid) &
172         ,stat=stat)
173    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
174    allocate(areaeast(0:numxgrid-1,0:numygrid-1,numzgrid) &
175         ,stat=stat)
176    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
177    allocate(areanorth(0:numxgrid-1,0:numygrid-1,numzgrid) &
178         ,stat=stat)
179    if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
180  return
181
182
183999   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "OUTGRID"    #### '
184  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
185  write(*,*) ' #### xxx/flexpart/options                    #### '
186  stop
187
188end subroutine readoutgrid
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG