1 | !*********************************************************************** |
---|
2 | !* Copyright 2012,2013 * |
---|
3 | !* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, * |
---|
4 | !* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,* |
---|
5 | !* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, * |
---|
6 | !* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, * |
---|
7 | !* * |
---|
8 | !* This file is part of FLEXPART WRF * |
---|
9 | !* * |
---|
10 | !* FLEXPART is free software: you can redistribute it and/or modify * |
---|
11 | !* it under the terms of the GNU General Public License as published by* |
---|
12 | !* the Free Software Foundation, either version 3 of the License, or * |
---|
13 | !* (at your option) any later version. * |
---|
14 | !* * |
---|
15 | !* FLEXPART is distributed in the hope that it will be useful, * |
---|
16 | !* but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
17 | !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
18 | !* GNU General Public License for more details. * |
---|
19 | !* * |
---|
20 | !* You should have received a copy of the GNU General Public License * |
---|
21 | !* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
22 | !*********************************************************************** |
---|
23 | |
---|
24 | subroutine openreceptors |
---|
25 | !******************************************************************************* |
---|
26 | ! * |
---|
27 | ! Note: This is the FLEXPART_WRF version of subroutine openreceptors. * |
---|
28 | ! * |
---|
29 | ! This routine opens the receptor output files and writes out the receptor * |
---|
30 | ! names and the receptor locations. The receptor output files are not closed, * |
---|
31 | ! but kept open throughout the simulation. Concentrations are continuously * |
---|
32 | ! dumped to these files. * |
---|
33 | ! * |
---|
34 | ! Author: A. Stohl * |
---|
35 | ! * |
---|
36 | ! 7 August 2002 * |
---|
37 | ! * |
---|
38 | ! Dec 2005, J. Fast - Output files can be either binary or ascii. * |
---|
39 | ! Write iomode_xycoord to output files. * |
---|
40 | ! Receptor positions can be lat-lon or grid-meters. * |
---|
41 | ! Dec 2005, R. Easter - changed names of "*lon0*" & "*lat0*" variables * |
---|
42 | ! * |
---|
43 | !******************************************************************************* |
---|
44 | ! * |
---|
45 | ! Variables: * |
---|
46 | ! numreceptor actual number of receptor points specified * |
---|
47 | ! receptornames names of the receptor points * |
---|
48 | ! xreceptor,yreceptor coordinates of the receptor points * |
---|
49 | ! * |
---|
50 | !******************************************************************************* |
---|
51 | |
---|
52 | use par_mod |
---|
53 | use com_mod |
---|
54 | |
---|
55 | implicit none |
---|
56 | |
---|
57 | integer :: j |
---|
58 | real :: xtmp(maxreceptor),ytmp(maxreceptor) |
---|
59 | real :: xtmpb,ytmpb |
---|
60 | |
---|
61 | |
---|
62 | ! Open output file for receptor points and write out a short header |
---|
63 | ! containing receptor names and locations |
---|
64 | !****************************************************************** |
---|
65 | |
---|
66 | if (numreceptor.ge.1) then ! do it only if receptors are specified |
---|
67 | |
---|
68 | do j = 1, numreceptor |
---|
69 | xtmp(j) = xreceptor(j)*dx + xmet0 |
---|
70 | ytmp(j) = yreceptor(j)*dy + ymet0 |
---|
71 | if (outgrid_option .eq. 0) then |
---|
72 | xtmpb = xtmp(j) |
---|
73 | ytmpb = ytmp(j) |
---|
74 | call xymeter_to_ll_wrf( xtmpb, ytmpb, xtmp(j), ytmp(j) ) |
---|
75 | endif |
---|
76 | enddo |
---|
77 | |
---|
78 | ! Concentration output |
---|
79 | !********************* |
---|
80 | |
---|
81 | if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then |
---|
82 | if (iouttype.eq.0) then |
---|
83 | open(unitoutrecept,file=path(1)(1:length(1))//'receptor_conc', & |
---|
84 | form='unformatted',err=997) |
---|
85 | write(unitoutrecept) (receptorname(j),j=1,numreceptor) |
---|
86 | write(unitoutrecept) (xtmp(j),ytmp(j),j=1,numreceptor), & |
---|
87 | outgrid_option |
---|
88 | endif |
---|
89 | if (iouttype.eq.1) then |
---|
90 | open(unitoutrecept,file=path(1)(1:length(1))//'receptor_conc', & |
---|
91 | form='formatted',err=997) |
---|
92 | ! do j = 1, numreceptor |
---|
93 | ! write(unitoutrecept,*) receptorname(j) |
---|
94 | ! enddo |
---|
95 | ! write(unitoutrecept,*) (xtmp(j),ytmp(j),j=1,numreceptor), |
---|
96 | do j=1,numreceptor |
---|
97 | write(unitoutrecept,1001) receptorname(j), & |
---|
98 | xtmp(j),ytmp(j),outgrid_option |
---|
99 | enddo |
---|
100 | 1001 format(a16,f10.4,f10.4,i5) |
---|
101 | endif |
---|
102 | endif |
---|
103 | |
---|
104 | ! Mixing ratio output |
---|
105 | !******************** |
---|
106 | |
---|
107 | if ((iout.eq.2).or.(iout.eq.3)) then |
---|
108 | if (iouttype.eq.0) then |
---|
109 | open(unitoutreceptppt,file=path(1)(1:length(1))//'receptor_pptv', & |
---|
110 | form='unformatted',err=998) |
---|
111 | write(unitoutreceptppt) (receptorname(j),j=1,numreceptor) |
---|
112 | write(unitoutreceptppt) (xtmp(j),ytmp(j),j=1,numreceptor), & |
---|
113 | outgrid_option |
---|
114 | endif |
---|
115 | if (iouttype.eq.1) then |
---|
116 | open(unitoutreceptppt,file=path(1)(1:length(1))//'receptor_pptv', & |
---|
117 | form='formatted',err=998) |
---|
118 | ! do j = 1, numreceptor |
---|
119 | ! write(unitoutreceptppt,*) receptorname(j) |
---|
120 | ! enddo |
---|
121 | ! write(unitoutreceptppt,*) (xtmp(j),ytmp(j),j=1,numreceptor), |
---|
122 | do j=1,numreceptor |
---|
123 | write(unitoutreceptppt,1001) receptorname(j), & |
---|
124 | xtmp(j),ytmp(j),outgrid_option |
---|
125 | enddo |
---|
126 | endif |
---|
127 | endif |
---|
128 | endif |
---|
129 | |
---|
130 | return |
---|
131 | |
---|
132 | |
---|
133 | 997 write(*,*) ' #### FLEXPART MODEL ERROR! THE FILE #### ' |
---|
134 | write(*,*) ' #### receptor_conc #### ' |
---|
135 | write(*,*) ' #### CANNOT BE OPENED. #### ' |
---|
136 | stop |
---|
137 | |
---|
138 | 998 write(*,*) ' #### FLEXPART MODEL ERROR! THE FILE #### ' |
---|
139 | write(*,*) ' #### receptor_pptv #### ' |
---|
140 | write(*,*) ' #### CANNOT BE OPENED. #### ' |
---|
141 | stop |
---|
142 | |
---|
143 | end subroutine openreceptors |
---|
144 | |
---|