[8] | 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 readreceptors |
---|
| 23 | |
---|
| 24 | !***************************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! This routine reads the user specifications for the receptor points. * |
---|
| 27 | ! * |
---|
| 28 | ! Author: A. Stohl * |
---|
| 29 | ! 1 August 1996 * |
---|
[10] | 30 | ! HSO, 14 August 2013 |
---|
| 31 | ! Added optional namelist input |
---|
[8] | 32 | ! * |
---|
| 33 | !***************************************************************************** |
---|
| 34 | ! * |
---|
| 35 | ! Variables: * |
---|
| 36 | ! receptorarea(maxreceptor) area of dx*dy at location of receptor * |
---|
| 37 | ! receptorname(maxreceptor) names of receptors * |
---|
| 38 | ! xreceptor,yreceptor coordinates of receptor points * |
---|
| 39 | ! * |
---|
| 40 | ! Constants: * |
---|
| 41 | ! unitreceptor unit connected to file RECEPTORS * |
---|
| 42 | ! * |
---|
| 43 | !***************************************************************************** |
---|
| 44 | |
---|
| 45 | use par_mod |
---|
| 46 | use com_mod |
---|
| 47 | |
---|
| 48 | implicit none |
---|
| 49 | |
---|
| 50 | integer :: j |
---|
| 51 | real :: x,y,xm,ym |
---|
| 52 | character(len=16) :: receptor |
---|
| 53 | |
---|
[10] | 54 | integer :: readerror |
---|
| 55 | real :: lon,lat ! for namelist input, lon/lat are used instead of x,y |
---|
[8] | 56 | |
---|
[10] | 57 | ! declare namelist |
---|
| 58 | namelist /receptors/ & |
---|
| 59 | receptor, lon, lat |
---|
| 60 | |
---|
| 61 | lon=-999.9 |
---|
| 62 | |
---|
[8] | 63 | ! For backward runs, do not allow receptor output. Thus, set number of receptors to zero |
---|
| 64 | !***************************************************************************** |
---|
| 65 | |
---|
| 66 | if (ldirect.lt.0) then |
---|
| 67 | numreceptor=0 |
---|
| 68 | return |
---|
| 69 | endif |
---|
| 70 | |
---|
| 71 | |
---|
| 72 | ! Open the RECEPTORS file and read output grid specifications |
---|
| 73 | !************************************************************ |
---|
| 74 | |
---|
| 75 | open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS', & |
---|
| 76 | status='old',err=999) |
---|
| 77 | |
---|
[10] | 78 | ! try namelist input |
---|
| 79 | read(unitreceptor,receptors,iostat=readerror) |
---|
[8] | 80 | |
---|
[10] | 81 | rewind(unitreceptor) |
---|
| 82 | if ((lon.lt.-900).or.(readerror.ne.0)) then |
---|
[8] | 83 | |
---|
[10] | 84 | call skplin(5,unitreceptor) |
---|
[8] | 85 | |
---|
[10] | 86 | ! Read the names and coordinates of the receptors |
---|
| 87 | !************************************************ |
---|
| 88 | |
---|
| 89 | j=0 |
---|
| 90 | 100 j=j+1 |
---|
[8] | 91 | read(unitreceptor,*,end=99) |
---|
| 92 | read(unitreceptor,*,end=99) |
---|
| 93 | read(unitreceptor,*,end=99) |
---|
| 94 | read(unitreceptor,'(4x,a16)',end=99) receptor |
---|
| 95 | call skplin(3,unitreceptor) |
---|
| 96 | read(unitreceptor,'(4x,f11.4)',end=99) x |
---|
| 97 | call skplin(3,unitreceptor) |
---|
| 98 | read(unitreceptor,'(4x,f11.4)',end=99) y |
---|
| 99 | if ((x.eq.0.).and.(y.eq.0.).and. & |
---|
| 100 | (receptor.eq.' ')) then |
---|
| 101 | j=j-1 |
---|
| 102 | goto 100 |
---|
| 103 | endif |
---|
| 104 | if (j.gt.maxreceptor) then |
---|
[10] | 105 | write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' |
---|
| 106 | write(*,*) ' #### POINTS ARE GIVEN. #### ' |
---|
| 107 | write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' |
---|
| 108 | write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' |
---|
[8] | 109 | endif |
---|
| 110 | receptorname(j)=receptor |
---|
| 111 | xreceptor(j)=(x-xlon0)/dx ! transform to grid coordinates |
---|
| 112 | yreceptor(j)=(y-ylat0)/dy |
---|
| 113 | xm=r_earth*cos(y*pi/180.)*dx/180.*pi |
---|
| 114 | ym=r_earth*dy/180.*pi |
---|
| 115 | receptorarea(j)=xm*ym |
---|
| 116 | goto 100 |
---|
| 117 | |
---|
[10] | 118 | 99 numreceptor=j-1 |
---|
| 119 | |
---|
| 120 | else ! continue with namelist input |
---|
| 121 | |
---|
| 122 | j=0 |
---|
| 123 | do while (readerror.eq.0) |
---|
| 124 | j=j+1 |
---|
| 125 | lon=-999.9 |
---|
| 126 | read(unitreceptor,receptors,iostat=readerror) |
---|
| 127 | if ((lon.lt.-900).or.(readerror.ne.0)) then |
---|
| 128 | readerror=1 |
---|
| 129 | else |
---|
| 130 | if (j.gt.maxreceptor) then |
---|
| 131 | write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' |
---|
| 132 | write(*,*) ' #### POINTS ARE GIVEN. #### ' |
---|
| 133 | write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' |
---|
| 134 | write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' |
---|
| 135 | endif |
---|
| 136 | receptorname(j)=receptor |
---|
| 137 | xreceptor(j)=(lon-xlon0)/dx ! transform to grid coordinates |
---|
| 138 | yreceptor(j)=(lat-ylat0)/dy |
---|
| 139 | xm=r_earth*cos(lat*pi/180.)*dx/180.*pi |
---|
| 140 | ym=r_earth*dy/180.*pi |
---|
| 141 | receptorarea(j)=xm*ym |
---|
| 142 | endif |
---|
| 143 | end do |
---|
| 144 | numreceptor=j-1 |
---|
| 145 | |
---|
| 146 | endif |
---|
| 147 | |
---|
| 148 | |
---|
[8] | 149 | close(unitreceptor) |
---|
| 150 | return |
---|
| 151 | |
---|
| 152 | |
---|
[10] | 153 | 999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' |
---|
[8] | 154 | write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' |
---|
| 155 | write(*,'(a)') path(1)(1:length(1)) |
---|
| 156 | stop |
---|
| 157 | |
---|
| 158 | end subroutine readreceptors |
---|