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 writeheader |
---|
25 | |
---|
26 | !***************************************************************************** |
---|
27 | ! * |
---|
28 | ! This routine produces a file header containing basic information on the * |
---|
29 | ! settings of the FLEXPART run. * |
---|
30 | ! The header file is essential and must be read in by any postprocessing * |
---|
31 | ! program before reading in the output data. * |
---|
32 | ! * |
---|
33 | ! Author: A. Stohl * |
---|
34 | ! * |
---|
35 | ! 7 August 2002 * |
---|
36 | ! * |
---|
37 | !***************************************************************************** |
---|
38 | ! * |
---|
39 | ! Variables: * |
---|
40 | ! * |
---|
41 | ! xlon longitude * |
---|
42 | ! xl model x coordinate * |
---|
43 | ! ylat latitude * |
---|
44 | ! yl model y coordinate * |
---|
45 | ! * |
---|
46 | !***************************************************************************** |
---|
47 | |
---|
48 | use point_mod |
---|
49 | use outg_mod |
---|
50 | use par_mod |
---|
51 | use com_mod |
---|
52 | |
---|
53 | implicit none |
---|
54 | |
---|
55 | integer :: jjjjmmdd,ihmmss,i,ix,jy,j |
---|
56 | real :: xp1,yp1,xp2,yp2 |
---|
57 | real :: xsw,xne,ysw,yne,tmpx,tmpy,tmplon,tmplat,xl2,yl2 |
---|
58 | |
---|
59 | |
---|
60 | !************************ |
---|
61 | ! Open header output file |
---|
62 | !************************ |
---|
63 | |
---|
64 | open(unitheader,file=path(1)(1:length(1))//'header', & |
---|
65 | form='unformatted',err=998) |
---|
66 | |
---|
67 | |
---|
68 | ! Write the header information |
---|
69 | !***************************** |
---|
70 | |
---|
71 | if (ldirect.eq.1) then |
---|
72 | ! write(unitheader) ibdate,ibtime,'FLEXWRF V2.1' |
---|
73 | if (outgrid_option .eq. 1) then |
---|
74 | write(unitheader) ibdate,ibtime,'FLEXWRF lalo ' |
---|
75 | else |
---|
76 | if (map_proj_id.eq.1) write(unitheader) ibdate,ibtime,'FLEXWRF lamb ' |
---|
77 | if (map_proj_id.eq.2) write(unitheader) ibdate,ibtime,'FLEXWRF ster ' |
---|
78 | if (map_proj_id.eq.3) write(unitheader) ibdate,ibtime,'FLEXWRF merc ' |
---|
79 | if (map_proj_id.eq.4) write(unitheader) ibdate,ibtime,'FLEXWRF glob ' |
---|
80 | endif |
---|
81 | else |
---|
82 | if (outgrid_option .eq. 1) then |
---|
83 | write(unitheader) iedate,ietime,'FLEXWRF lalo ' |
---|
84 | else |
---|
85 | if (map_proj_id.eq.1) write(unitheader) iedate,ietime,'FLEXWRF lamb ' |
---|
86 | if (map_proj_id.eq.2) write(unitheader) iedate,ietime,'FLEXWRF ster ' |
---|
87 | if (map_proj_id.eq.3) write(unitheader) iedate,ietime,'FLEXWRF merc ' |
---|
88 | if (map_proj_id.eq.4) write(unitheader) iedate,ietime,'FLEXWRF glob ' |
---|
89 | endif |
---|
90 | endif |
---|
91 | |
---|
92 | ! Write info on output interval, averaging time, sampling time |
---|
93 | !************************************************************* |
---|
94 | |
---|
95 | write(unitheader) loutstep,loutaver,loutsample |
---|
96 | |
---|
97 | ! Write information on output grid setup |
---|
98 | !*************************************** |
---|
99 | |
---|
100 | if (outgrid_option .eq. 1) then |
---|
101 | write(unitheader) outlon0,outlat0,numxgrid,numygrid, & |
---|
102 | dxoutl,dyoutl |
---|
103 | else |
---|
104 | write(unitheader) outlon0,outlat0,numxgrid,numygrid, & |
---|
105 | dxout,dyout |
---|
106 | endif |
---|
107 | write(unitheader) numzgrid,(outheight(i),i=1,numzgrid) |
---|
108 | |
---|
109 | call caldate(bdate,jjjjmmdd,ihmmss) |
---|
110 | write(unitheader) jjjjmmdd,ihmmss |
---|
111 | |
---|
112 | ! Write number of species, and name for each species (+extra name for depositions) |
---|
113 | ! Indicate the dimension of the fields (i.e., 1 for deposition fields, numzgrid for |
---|
114 | ! concentration fields |
---|
115 | !***************************************************************************** |
---|
116 | |
---|
117 | write(unitheader) 3*nspec,maxpointspec_act |
---|
118 | do i=1,nspec |
---|
119 | write(unitheader) 1,'WD_'//species(i)(1:7) |
---|
120 | write(unitheader) 1,'DD_'//species(i)(1:7) |
---|
121 | write(unitheader) numzgrid,species(i) |
---|
122 | end do |
---|
123 | |
---|
124 | ! Write information on release points: total number, then for each point: |
---|
125 | ! start, end, coordinates, # of particles, name, mass |
---|
126 | !************************************************************************ |
---|
127 | |
---|
128 | write(unitheader) numpoint |
---|
129 | do i=1,numpoint |
---|
130 | write(unitheader) ireleasestart(i),ireleaseend(i),kindz(i) |
---|
131 | xp1=xpoint1(i)*dx+xlon0 |
---|
132 | yp1=ypoint1(i)*dy+ylat0 |
---|
133 | xp2=xpoint2(i)*dx+xlon0 |
---|
134 | yp2=ypoint2(i)*dy+ylat0 |
---|
135 | write(unitheader) xp1,yp1,xp2,yp2,zpoint1(i),zpoint2(i) |
---|
136 | write(unitheader) npart(i),1 |
---|
137 | if (numpoint.le.2000) then |
---|
138 | write(unitheader) compoint(i) |
---|
139 | else |
---|
140 | write(unitheader) compoint(2001) |
---|
141 | endif |
---|
142 | do j=1,nspec |
---|
143 | write(unitheader) xmass(i,j) |
---|
144 | write(unitheader) xmass(i,j) |
---|
145 | write(unitheader) xmass(i,j) |
---|
146 | end do |
---|
147 | end do |
---|
148 | |
---|
149 | ! Write information on some model switches |
---|
150 | !***************************************** |
---|
151 | |
---|
152 | write(unitheader) method,lsubgrid,lconvection, & |
---|
153 | ind_source,ind_receptor |
---|
154 | |
---|
155 | ! Write age class information |
---|
156 | !**************************** |
---|
157 | |
---|
158 | write(unitheader) nageclass,(lage(i),i=1,nageclass) |
---|
159 | |
---|
160 | |
---|
161 | ! Write topography to output file |
---|
162 | !******************************** |
---|
163 | |
---|
164 | do ix=0,numxgrid-1 |
---|
165 | write(unitheader) (oroout(ix,jy),jy=0,numygrid-1) |
---|
166 | end do |
---|
167 | close(unitheader) |
---|
168 | |
---|
169 | |
---|
170 | open(53,file=path(1)(1:length(1))//'latlon.txt',form='formatted') |
---|
171 | open(54,file=path(1)(1:length(1))//'latlon_corner.txt' & |
---|
172 | ,form='formatted') |
---|
173 | |
---|
174 | if (outgrid_option.eq.0) then ! irregular |
---|
175 | call ll_to_xymeter_wrf(outgrid_swlon,outgrid_swlat,xsw,ysw) |
---|
176 | call ll_to_xymeter_wrf(outgrid_nelon,outgrid_nelat,xne,yne) |
---|
177 | do jy=1,numygrid |
---|
178 | do ix=1,numxgrid |
---|
179 | tmpx=out_xm0+(float(ix)-0.5)*dxout |
---|
180 | tmpy=out_ym0+(float(jy)-0.5)*dyout |
---|
181 | call xymeter_to_ll_wrf(tmpx,tmpy,tmplon,tmplat) |
---|
182 | write(53,*) tmplon,tmplat |
---|
183 | tmpx=out_xm0+(float(ix)-1.)*dxout |
---|
184 | tmpy=out_ym0+(float(jy)-1.)*dyout |
---|
185 | call xymeter_to_ll_wrf_out(tmpx,tmpy,tmplon,tmplat) |
---|
186 | write(54,*) tmplon,tmplat |
---|
187 | enddo |
---|
188 | enddo |
---|
189 | else ! regular |
---|
190 | call ll_to_xymeter_wrf(outgrid_swlon,outgrid_swlat,xsw,ysw) |
---|
191 | call ll_to_xymeter_wrf(outgrid_nelon,outgrid_nelat,xne,yne) |
---|
192 | do jy=1,numygrid |
---|
193 | do ix=1,numxgrid |
---|
194 | tmpx=xsw+(xne-xsw)*float(ix-1)/float(numxgrid-1) |
---|
195 | tmpy=ysw+(yne-ysw)*float(jy-1)/float(numygrid-1) |
---|
196 | call xymeter_to_ll_wrf(tmpx,tmpy,tmplon,tmplat) |
---|
197 | xl2=outlon0+(float(ix)-0.5)*dxoutl !long |
---|
198 | yl2=outlat0+(float(jy)-0.5)*dyoutl !lat |
---|
199 | write(53,*) xl2,yl2 |
---|
200 | xl2=outlon0+float(ix-1)*dxoutl !long |
---|
201 | yl2=outlat0+float(jy-1)*dyoutl !lat |
---|
202 | write(54,*) xl2,yl2 |
---|
203 | enddo |
---|
204 | enddo |
---|
205 | endif |
---|
206 | |
---|
207 | |
---|
208 | close(53) |
---|
209 | close(54) |
---|
210 | |
---|
211 | return |
---|
212 | |
---|
213 | 998 write(*,*) ' #### FLEXPART MODEL ERROR! THE FILE #### ' |
---|
214 | write(*,*) ' #### '//path(1)(1:length(1))//'header'//' #### ' |
---|
215 | write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS #### ' |
---|
216 | write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### ' |
---|
217 | write(*,*) ' #### THE PROGRAM AGAIN. #### ' |
---|
218 | stop |
---|
219 | |
---|
220 | end subroutine writeheader |
---|