CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
readoutgrid_nest.f90
Go to the documentation of this file.
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 readoutgrid_nest
23 
24  !*****************************************************************************
25  ! *
26  ! This routine reads the user specifications for the output nest. *
27  ! *
28  ! Author: A. Stohl *
29  ! *
30  ! 4 June 1996 *
31  ! *
32  !*****************************************************************************
33  ! *
34  ! Variables: *
35  ! dxoutn,dyoutn grid distances of output nest *
36  ! numxgridn,numygridn,numzgrid nest dimensions *
37  ! outlon0n,outlat0n lower left corner of nest *
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 :: stat
52  real :: 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_NEST', &
61  status='old',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)') outlon0n
72  call skplin(3,unitoutgrid)
73  read(unitoutgrid,'(4x,f11.4)') outlat0n
74  call skplin(3,unitoutgrid)
75  read(unitoutgrid,'(4x,i5)') numxgridn
76  call skplin(3,unitoutgrid)
77  read(unitoutgrid,'(4x,i5)') numygridn
78  call skplin(3,unitoutgrid)
79  read(unitoutgrid,'(4x,f12.5)') dxoutn
80  call skplin(3,unitoutgrid)
81  read(unitoutgrid,'(4x,f12.5)') dyoutn
82 
83 
84  allocate(orooutn(0:numxgridn-1,0:numygridn-1) &
85  ,stat=stat)
86  if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
87  allocate(arean(0:numxgridn-1,0:numygridn-1) &
88  ,stat=stat)
89  if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
90  allocate(volumen(0:numxgridn-1,0:numygridn-1,numzgrid) &
91  ,stat=stat)
92  if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
93 
94  ! Check validity of output grid (shall be within model domain)
95  !*************************************************************
96 
97  xr=outlon0n+real(numxgridn)*dxoutn
98  yr=outlat0n+real(numygridn)*dyoutn
99  xr1=xlon0+real(nxmin1)*dx
100  yr1=ylat0+real(nymin1)*dy
101  if ((outlon0n+eps.lt.xlon0).or.(outlat0n+eps.lt.ylat0) &
102  .or.(xr.gt.xr1+eps).or.(yr.gt.yr1+eps)) then
103  write(*,*) ' #### FLEXPART MODEL ERROR! PART OF OUTPUT ####'
104  write(*,*) ' #### NEST IS OUTSIDE MODEL DOMAIN. CHANGE ####'
105  write(*,*) ' #### FILE OUTGRID IN DIRECTORY ####'
106  write(*,'(a)') path(1)(1:length(1))
107  stop
108  endif
109 
110  xoutshiftn=xlon0-outlon0n
111  youtshiftn=ylat0-outlat0n
112  close(unitoutgrid)
113  return
114 
115 
116 999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE OUTGRID_NEST #### '
117  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### '
118  write(*,*) ' #### xxx/flexpart/options #### '
119  stop
120 
121 end subroutine readoutgrid_nest
subroutine skplin(nlines, iunit)
Definition: skplin.f90:22
subroutine readoutgrid_nest