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 * |
---|
30 | ! HSO, 14 August 2013 |
---|
31 | ! Added optional namelist input |
---|
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 | |
---|
54 | integer :: readerror |
---|
55 | real :: lon,lat ! for namelist input, lon/lat are used instead of x,y |
---|
56 | integer,parameter :: unitreceptorout=2 |
---|
57 | |
---|
58 | ! declare namelist |
---|
59 | namelist /receptors/ & |
---|
60 | receptor, lon, lat |
---|
61 | |
---|
62 | lon=-999.9 |
---|
63 | |
---|
64 | ! For backward runs, do not allow receptor output. Thus, set number of receptors to zero |
---|
65 | !***************************************************************************** |
---|
66 | |
---|
67 | if (ldirect.lt.0) then |
---|
68 | numreceptor=0 |
---|
69 | return |
---|
70 | endif |
---|
71 | |
---|
72 | |
---|
73 | ! Open the RECEPTORS file and read output grid specifications |
---|
74 | !************************************************************ |
---|
75 | |
---|
76 | open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',form='formatted',status='old',err=999) |
---|
77 | |
---|
78 | ! try namelist input |
---|
79 | read(unitreceptor,receptors,iostat=readerror) |
---|
80 | |
---|
81 | ! prepare namelist output if requested |
---|
82 | if (nmlout.and.lroot) then |
---|
83 | open(unitreceptorout,file=path(2)(1:length(2))//'RECEPTORS.namelist',& |
---|
84 | &access='append',status='replace',err=1000) |
---|
85 | endif |
---|
86 | |
---|
87 | if ((lon.lt.-900).or.(readerror.ne.0)) then |
---|
88 | |
---|
89 | close(unitreceptor) |
---|
90 | open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',status='old',err=999) |
---|
91 | call skplin(5,unitreceptor) |
---|
92 | |
---|
93 | ! Read the names and coordinates of the receptors |
---|
94 | !************************************************ |
---|
95 | |
---|
96 | j=0 |
---|
97 | 100 j=j+1 |
---|
98 | read(unitreceptor,*,end=99) |
---|
99 | read(unitreceptor,*,end=99) |
---|
100 | read(unitreceptor,*,end=99) |
---|
101 | read(unitreceptor,'(4x,a16)',end=99) receptor |
---|
102 | call skplin(3,unitreceptor) |
---|
103 | read(unitreceptor,'(4x,f11.4)',end=99) x |
---|
104 | call skplin(3,unitreceptor) |
---|
105 | read(unitreceptor,'(4x,f11.4)',end=99) y |
---|
106 | if ((x.eq.0.).and.(y.eq.0.).and. & |
---|
107 | (receptor.eq.' ')) then |
---|
108 | j=j-1 |
---|
109 | goto 100 |
---|
110 | endif |
---|
111 | if (j.gt.maxreceptor) then |
---|
112 | write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' |
---|
113 | write(*,*) ' #### POINTS ARE GIVEN. #### ' |
---|
114 | write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' |
---|
115 | write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' |
---|
116 | endif |
---|
117 | receptorname(j)=receptor |
---|
118 | xreceptor(j)=(x-xlon0)/dx ! transform to grid coordinates |
---|
119 | yreceptor(j)=(y-ylat0)/dy |
---|
120 | xm=r_earth*cos(y*pi/180.)*dx/180.*pi |
---|
121 | ym=r_earth*dy/180.*pi |
---|
122 | receptorarea(j)=xm*ym |
---|
123 | |
---|
124 | ! write receptors file in namelist format to output directory if requested |
---|
125 | if (nmlout.and.lroot) then |
---|
126 | lon=x |
---|
127 | lat=y |
---|
128 | write(unitreceptorout,nml=receptors) |
---|
129 | endif |
---|
130 | |
---|
131 | goto 100 |
---|
132 | |
---|
133 | 99 numreceptor=j-1 |
---|
134 | |
---|
135 | else ! continue with namelist input |
---|
136 | |
---|
137 | j=0 |
---|
138 | do while (readerror.eq.0) |
---|
139 | j=j+1 |
---|
140 | lon=-999.9 |
---|
141 | read(unitreceptor,receptors,iostat=readerror) |
---|
142 | if ((lon.lt.-900).or.(readerror.ne.0)) then |
---|
143 | readerror=1 |
---|
144 | else |
---|
145 | if (j.gt.maxreceptor) then |
---|
146 | write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' |
---|
147 | write(*,*) ' #### POINTS ARE GIVEN. #### ' |
---|
148 | write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' |
---|
149 | write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' |
---|
150 | endif |
---|
151 | receptorname(j)=receptor |
---|
152 | xreceptor(j)=(lon-xlon0)/dx ! transform to grid coordinates |
---|
153 | yreceptor(j)=(lat-ylat0)/dy |
---|
154 | xm=r_earth*cos(lat*pi/180.)*dx/180.*pi |
---|
155 | ym=r_earth*dy/180.*pi |
---|
156 | receptorarea(j)=xm*ym |
---|
157 | endif |
---|
158 | |
---|
159 | ! write receptors file in namelist format to output directory if requested |
---|
160 | if (nmlout.and.lroot) then |
---|
161 | write(unitreceptorout,nml=receptors) |
---|
162 | endif |
---|
163 | |
---|
164 | end do |
---|
165 | numreceptor=j-1 |
---|
166 | close(unitreceptor) |
---|
167 | |
---|
168 | endif |
---|
169 | |
---|
170 | if (nmlout.and.lroot) then |
---|
171 | close(unitreceptorout) |
---|
172 | endif |
---|
173 | |
---|
174 | return |
---|
175 | |
---|
176 | |
---|
177 | 999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' |
---|
178 | write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' |
---|
179 | write(*,'(a)') path(1)(1:length(1)) |
---|
180 | stop |
---|
181 | |
---|
182 | 1000 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' |
---|
183 | write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' |
---|
184 | write(*,'(a)') path(2)(1:length(2)) |
---|
185 | stop |
---|
186 | |
---|
187 | end subroutine readreceptors |
---|