1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121
|
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann *
! *
! This file is part of FLEXPART. *
! *
! FLEXPART is free software: you can redistribute it and/or modify *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or *
! (at your option) any later version. *
! *
! FLEXPART is distributed in the hope that it will be useful, *
! but WITHOUT ANY WARRANTY; without even the implied warranty of *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
! GNU General Public License for more details. *
! *
! You should have received a copy of the GNU General Public License *
! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
!**********************************************************************
subroutine readoutgrid_nest
!*****************************************************************************
! *
! This routine reads the user specifications for the output nest. *
! *
! Author: A. Stohl *
! *
! 4 June 1996 *
! *
!*****************************************************************************
! *
! Variables: *
! dxoutn,dyoutn grid distances of output nest *
! numxgridn,numygridn,numzgrid nest dimensions *
! outlon0n,outlat0n lower left corner of nest *
! outheight(maxzgrid) height levels of output grid [m] *
! *
! Constants: *
! unitoutgrid unit connected to file OUTGRID *
! *
!*****************************************************************************
use outg_mod
use par_mod
use com_mod
implicit none
integer :: stat
real :: xr,xr1,yr,yr1
real,parameter :: eps=1.e-4
! Open the OUTGRID file and read output grid specifications
!**********************************************************
open(unitoutgrid,file=path(1)(1:length(1))//'OUTGRID_NEST', &
status='old',err=999)
call skplin(5,unitoutgrid)
! 1. Read horizontal grid specifications
!****************************************
call skplin(3,unitoutgrid)
read(unitoutgrid,'(4x,f11.4)') outlon0n
call skplin(3,unitoutgrid)
read(unitoutgrid,'(4x,f11.4)') outlat0n
call skplin(3,unitoutgrid)
read(unitoutgrid,'(4x,i5)') numxgridn
call skplin(3,unitoutgrid)
read(unitoutgrid,'(4x,i5)') numygridn
call skplin(3,unitoutgrid)
read(unitoutgrid,'(4x,f12.5)') dxoutn
call skplin(3,unitoutgrid)
read(unitoutgrid,'(4x,f12.5)') dyoutn
allocate(orooutn(0:numxgridn-1,0:numygridn-1) &
,stat=stat)
if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
allocate(arean(0:numxgridn-1,0:numygridn-1) &
,stat=stat)
if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
allocate(volumen(0:numxgridn-1,0:numygridn-1,numzgrid) &
,stat=stat)
if (stat.ne.0) write(*,*)'ERROR: could not allocate outh'
! Check validity of output grid (shall be within model domain)
!*************************************************************
xr=outlon0n+real(numxgridn)*dxoutn
yr=outlat0n+real(numygridn)*dyoutn
xr1=xlon0+real(nxmin1)*dx
yr1=ylat0+real(nymin1)*dy
if ((outlon0n+eps.lt.xlon0).or.(outlat0n+eps.lt.ylat0) &
.or.(xr.gt.xr1+eps).or.(yr.gt.yr1+eps)) then
write(*,*) ' #### FLEXPART MODEL ERROR! PART OF OUTPUT ####'
write(*,*) ' #### NEST IS OUTSIDE MODEL DOMAIN. CHANGE ####'
write(*,*) ' #### FILE OUTGRID IN DIRECTORY ####'
write(*,'(a)') path(1)(1:length(1))
stop
endif
xoutshiftn=xlon0-outlon0n
youtshiftn=ylat0-outlat0n
close(unitoutgrid)
return
999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE OUTGRID_NEST #### '
write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### '
write(*,*) ' #### xxx/flexpart/options #### '
stop
end subroutine readoutgrid_nest
|