1 | program streamfunction |
---|
2 | |
---|
3 | ! SL 12/2009: |
---|
4 | ! This program reads 4D (lon-lat-alt-time) fields directly from LMD outputs |
---|
5 | ! without regrid : histmth OR from files recast in log P coordinates (_P) |
---|
6 | ! |
---|
7 | ! it computes: |
---|
8 | ! dmass -- 4D -- mass of each cell |
---|
9 | ! psi -- 3D -- Stream function |
---|
10 | ! |
---|
11 | ! Minimal requirements and dependencies: |
---|
12 | ! The dataset must include the following data: |
---|
13 | ! - surface pressure and surface geopotential |
---|
14 | ! - meridional winds |
---|
15 | |
---|
16 | implicit none |
---|
17 | |
---|
18 | include "netcdf.inc" ! NetCDF definitions |
---|
19 | |
---|
20 | character (len=128) :: infile ! input file name (name_P.nc) |
---|
21 | character (len=128) :: outfile ! output file name |
---|
22 | |
---|
23 | character (len=64) :: text ! to store some text |
---|
24 | integer infid ! NetCDF input file ID |
---|
25 | integer outfid ! NetCDF output file ID |
---|
26 | integer lon_dimid,lat_dimid,alt_dimid,time_dimid ! NetCDF dimension IDs |
---|
27 | integer lon_varid,lat_varid,alt_varid,time_varid |
---|
28 | integer :: datashape1d ! shape of 1D datasets |
---|
29 | integer,dimension(2) :: datashape2d ! shape of 2D datasets |
---|
30 | integer,dimension(3) :: datashape3d ! shape of 3D datasets |
---|
31 | integer,dimension(4) :: datashape4d ! shape of 4D datasets |
---|
32 | |
---|
33 | real :: miss_val=9.99e+33 ! special "missing value" to specify missing data |
---|
34 | real,parameter :: miss_val_def=9.99e+33 ! default value for "missing value" |
---|
35 | real :: pi |
---|
36 | real,dimension(:),allocatable :: lon ! longitude |
---|
37 | integer lonlength ! # of grid points along longitude |
---|
38 | real,dimension(:),allocatable :: lat ! latitude |
---|
39 | real,dimension(:),allocatable :: latrad ! latitude in radian |
---|
40 | integer latlength ! # of grid points along latitude |
---|
41 | real,dimension(:),allocatable :: plev ! Pressure levels (Pa) |
---|
42 | integer altlength ! # of grid point along altitude (of input datasets) |
---|
43 | real,dimension(:),allocatable :: time ! time |
---|
44 | integer timelength ! # of points along time |
---|
45 | real,dimension(:,:,:),allocatable :: ps ! surface pressure |
---|
46 | real,dimension(:,:,:,:),allocatable :: vitv ! meridional wind (in m/s) |
---|
47 | real,dimension(:,:,:,:),allocatable :: za ! areoid levels (m) |
---|
48 | |
---|
49 | real,dimension(:,:,:,:),allocatable :: rayon ! distance to center (m) |
---|
50 | real,dimension(:,:,:,:),allocatable :: grav ! gravity field (m s-2) |
---|
51 | real,dimension(:,:,:,:),allocatable :: dmass ! mass in cell (kg) |
---|
52 | real,dimension(:,:,:,:),allocatable :: vm ! meridional mass flux |
---|
53 | real,dimension(:,:,:),allocatable :: psi ! stream function |
---|
54 | |
---|
55 | integer ierr,ierr1,ierr2 ! NetCDF routines return codes |
---|
56 | integer i,j,ilon,ilat,ilev,itim ! for loops |
---|
57 | logical :: lmdflag |
---|
58 | |
---|
59 | real :: deltalat,deltalon ! lat and lon intervals in radians |
---|
60 | |
---|
61 | include "planet.h" |
---|
62 | |
---|
63 | !=============================================================================== |
---|
64 | ! 1. Input parameters |
---|
65 | !=============================================================================== |
---|
66 | |
---|
67 | pi = 2.*asin(1.) |
---|
68 | |
---|
69 | write(*,*) "" |
---|
70 | write(*,*) "You are working on the atmosphere of ",planet |
---|
71 | |
---|
72 | !=============================================================================== |
---|
73 | ! 1.1 Input file |
---|
74 | !=============================================================================== |
---|
75 | |
---|
76 | write(*,*) "" |
---|
77 | write(*,*) " Program valid for Venus or Titan LMD output files" |
---|
78 | write(*,*) "Enter input file name:" |
---|
79 | |
---|
80 | read(*,'(a128)') infile |
---|
81 | write(*,*) "" |
---|
82 | |
---|
83 | ! open input file |
---|
84 | |
---|
85 | ierr = NF_OPEN(infile,NF_NOWRITE,infid) |
---|
86 | if (ierr.ne.NF_NOERR) then |
---|
87 | write(*,*) 'ERROR: Pb opening file ',trim(infile) |
---|
88 | stop "" |
---|
89 | endif |
---|
90 | |
---|
91 | !=============================================================================== |
---|
92 | ! 1.2 Get grids in lon,lat,alt(pressure),time |
---|
93 | !=============================================================================== |
---|
94 | |
---|
95 | call get_iddim(infid,lat_varid,latlength,lon_varid,lonlength,& |
---|
96 | alt_varid,altlength,time_varid,timelength,lmdflag ) |
---|
97 | |
---|
98 | allocate(lon(lonlength)) |
---|
99 | ierr=NF_GET_VAR_REAL(infid,lon_varid,lon) |
---|
100 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading longitude" |
---|
101 | |
---|
102 | allocate(lat(latlength)) |
---|
103 | ierr=NF_GET_VAR_REAL(infid,lat_varid,lat) |
---|
104 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading lat" |
---|
105 | |
---|
106 | allocate(latrad(latlength)) |
---|
107 | latrad = lat*pi/180. |
---|
108 | |
---|
109 | ! Lat, lon pressure intervals |
---|
110 | deltalat = abs(latrad(2)-latrad(1)) |
---|
111 | deltalon = abs(lon(2)-lon(1))*pi/180. |
---|
112 | |
---|
113 | allocate(plev(altlength)) |
---|
114 | ierr=NF_GET_VAR_REAL(infid,alt_varid,plev) |
---|
115 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading altitude (ie pressure levels)" |
---|
116 | |
---|
117 | allocate(time(timelength)) |
---|
118 | ierr=NF_GET_VAR_REAL(infid,time_varid,time) |
---|
119 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading time" |
---|
120 | |
---|
121 | ! Time axis IN PLANET DAYS |
---|
122 | |
---|
123 | time=time/localday |
---|
124 | |
---|
125 | !=============================================================================== |
---|
126 | ! 1.3 Get output file name |
---|
127 | !=============================================================================== |
---|
128 | write(*,*) "" |
---|
129 | !write(*,*) "Enter output file name" |
---|
130 | !read(*,*) outfile |
---|
131 | outfile=infile(1:len_trim(infile)-3)//"_PSI.nc" |
---|
132 | write(*,*) "Output file name is: "//trim(outfile) |
---|
133 | |
---|
134 | |
---|
135 | |
---|
136 | !=============================================================================== |
---|
137 | ! 2.1 Store needed fields |
---|
138 | !=============================================================================== |
---|
139 | |
---|
140 | !=============================================================================== |
---|
141 | ! 2.1.1 Surface pressure |
---|
142 | !=============================================================================== |
---|
143 | allocate(ps(lonlength,latlength,timelength)) |
---|
144 | |
---|
145 | text="ps" |
---|
146 | call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2) |
---|
147 | if (ierr1.ne.NF_NOERR) then |
---|
148 | write(*,*) " looking for psol instead... " |
---|
149 | text="psol" |
---|
150 | call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2) |
---|
151 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get psol ID" |
---|
152 | endif |
---|
153 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface pressure" |
---|
154 | |
---|
155 | !=============================================================================== |
---|
156 | ! 2.1.3 Winds |
---|
157 | !=============================================================================== |
---|
158 | allocate(vitv(lonlength,latlength,altlength,timelength)) |
---|
159 | |
---|
160 | ! meridional wind vitv (in m/s) |
---|
161 | text="vitv" |
---|
162 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitv,ierr1,ierr2) |
---|
163 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitv ID" |
---|
164 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading meridional wind" |
---|
165 | |
---|
166 | !=============================================================================== |
---|
167 | ! 2.1.4 Altitude above areoide |
---|
168 | !=============================================================================== |
---|
169 | |
---|
170 | allocate(za(lonlength,latlength,altlength,timelength)) |
---|
171 | ! present only in _P regrided files |
---|
172 | ! For others, using g0*a0*a0/(g0*a0-geop)-a0 |
---|
173 | |
---|
174 | text="zareoid" |
---|
175 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,za,ierr1,ierr2) |
---|
176 | if (ierr1.ne.NF_NOERR) then |
---|
177 | write(*,*) " looking for geop instead... " |
---|
178 | text="geop" |
---|
179 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,za,ierr1,ierr2) |
---|
180 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get geop ID" |
---|
181 | do itim=1,timelength |
---|
182 | do ilon=1,lonlength |
---|
183 | do ilat=1,latlength |
---|
184 | do ilev=1,altlength |
---|
185 | if (za(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
186 | za(ilon,ilat,ilev,itim) = (g0*a0*a0/(g0*a0-za(ilon,ilat,ilev,itim))-a0)/1000. ! in km |
---|
187 | else |
---|
188 | za(ilon,ilat,ilev,itim) = miss_val |
---|
189 | endif |
---|
190 | enddo |
---|
191 | enddo |
---|
192 | enddo |
---|
193 | enddo |
---|
194 | endif |
---|
195 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zareoid/geop" |
---|
196 | |
---|
197 | !=============================================================================== |
---|
198 | ! 2.2 Computations |
---|
199 | !=============================================================================== |
---|
200 | |
---|
201 | !=============================================================================== |
---|
202 | ! 2.2.2 Mass in cells |
---|
203 | !=============================================================================== |
---|
204 | allocate(rayon(lonlength,latlength,altlength,timelength)) |
---|
205 | allocate(grav(lonlength,latlength,altlength,timelength)) |
---|
206 | allocate(dmass(lonlength,latlength,altlength,timelength)) |
---|
207 | |
---|
208 | do itim=1,timelength |
---|
209 | do ilon=1,lonlength |
---|
210 | do ilat=1,latlength |
---|
211 | do ilev=1,altlength |
---|
212 | ! Need to be consistent with GCM computations |
---|
213 | if (za(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
214 | rayon(ilon,ilat,ilev,itim) = a0 |
---|
215 | ! rayon(ilon,ilat,ilev,itim) = za(ilon,ilat,ilev,itim) + a0 |
---|
216 | grav(ilon,ilat,ilev,itim) = g0*a0*a0 & |
---|
217 | /(rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim)) |
---|
218 | else |
---|
219 | rayon(ilon,ilat,ilev,itim) = miss_val |
---|
220 | grav(ilon,ilat,ilev,itim) = miss_val |
---|
221 | endif |
---|
222 | enddo |
---|
223 | enddo |
---|
224 | enddo |
---|
225 | enddo ! timelength |
---|
226 | |
---|
227 | call cellmass(infid,latlength,lonlength,altlength,timelength, & |
---|
228 | lmdflag,deltalon,deltalat,latrad,plev,ps,grav,rayon, & |
---|
229 | dmass ) |
---|
230 | |
---|
231 | !=============================================================================== |
---|
232 | ! 2.2.8 Stream function |
---|
233 | !=============================================================================== |
---|
234 | allocate(vm(lonlength,latlength,altlength,timelength)) |
---|
235 | allocate(psi(latlength,altlength,timelength)) |
---|
236 | |
---|
237 | do itim=1,timelength |
---|
238 | |
---|
239 | do ilat=1,latlength |
---|
240 | do ilev=1,altlength |
---|
241 | do ilon=1,lonlength |
---|
242 | if (dmass(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
243 | vm(ilon,ilat,ilev,itim) = vitv(ilon,ilat,ilev,itim) & |
---|
244 | * dmass(ilon,ilat,ilev,itim) & |
---|
245 | / (rayon(ilon,ilat,ilev,itim)*deltalat) |
---|
246 | else |
---|
247 | vm(ilon,ilat,ilev,itim) = miss_val |
---|
248 | endif |
---|
249 | enddo |
---|
250 | enddo |
---|
251 | enddo |
---|
252 | |
---|
253 | |
---|
254 | do ilat=1,latlength |
---|
255 | psi(ilat,altlength,itim) = 0. |
---|
256 | do ilon=1,lonlength |
---|
257 | if (vm(ilon,ilat,altlength,itim).ne.miss_val) then |
---|
258 | psi(ilat,altlength,itim) = psi(ilat,altlength,itim) & |
---|
259 | + vm(ilon,ilat,altlength,itim) |
---|
260 | endif |
---|
261 | enddo |
---|
262 | do ilev=altlength-1,1,-1 |
---|
263 | psi(ilat,ilev,itim) = psi(ilat,ilev+1,itim) |
---|
264 | do ilon=1,lonlength |
---|
265 | if (vm(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
266 | psi(ilat,ilev,itim) = psi(ilat,ilev,itim) & |
---|
267 | + vm(ilon,ilat,ilev,itim) |
---|
268 | endif |
---|
269 | enddo |
---|
270 | enddo |
---|
271 | enddo |
---|
272 | |
---|
273 | enddo ! timelength |
---|
274 | |
---|
275 | print*,"End of computations" |
---|
276 | |
---|
277 | !=============================================================================== |
---|
278 | ! 3. Create output file |
---|
279 | !=============================================================================== |
---|
280 | |
---|
281 | ! Create output file |
---|
282 | ierr=NF_CREATE(outfile,NF_CLOBBER,outfid) |
---|
283 | if (ierr.ne.NF_NOERR) then |
---|
284 | write(*,*)"Error: could not create file ",outfile |
---|
285 | stop |
---|
286 | endif |
---|
287 | |
---|
288 | !=============================================================================== |
---|
289 | ! 3.1. Define and write dimensions |
---|
290 | !=============================================================================== |
---|
291 | |
---|
292 | call write_dim(outfid,lonlength,latlength,altlength,timelength, & |
---|
293 | lon,lat,plev,time,lon_dimid,lat_dimid,alt_dimid,time_dimid) |
---|
294 | |
---|
295 | !=============================================================================== |
---|
296 | ! 3.2. Define and write variables |
---|
297 | !=============================================================================== |
---|
298 | |
---|
299 | ! 1D Variables |
---|
300 | |
---|
301 | datashape1d=time_dimid |
---|
302 | |
---|
303 | ! 3D variables |
---|
304 | |
---|
305 | datashape3d(1)=lon_dimid |
---|
306 | datashape3d(2)=lat_dimid |
---|
307 | datashape3d(3)=time_dimid |
---|
308 | |
---|
309 | call write_var3d(outfid,datashape3d,lonlength,latlength,timelength,& |
---|
310 | "ps ", "Surface pressure ","Pa ",miss_val,& |
---|
311 | ps ) |
---|
312 | |
---|
313 | datashape3d(1)=lat_dimid |
---|
314 | datashape3d(2)=alt_dimid |
---|
315 | datashape3d(3)=time_dimid |
---|
316 | |
---|
317 | call write_var3d(outfid,datashape3d,lonlength,latlength,timelength,& |
---|
318 | "psi ", "Stream function ","kg/s ",miss_val,& |
---|
319 | psi ) |
---|
320 | |
---|
321 | ! 4D variables |
---|
322 | |
---|
323 | datashape4d(1)=lon_dimid |
---|
324 | datashape4d(2)=lat_dimid |
---|
325 | datashape4d(3)=alt_dimid |
---|
326 | datashape4d(4)=time_dimid |
---|
327 | |
---|
328 | call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
329 | "dmass ", "Mass ","kg ",miss_val,& |
---|
330 | dmass ) |
---|
331 | |
---|
332 | !!!! Close output file |
---|
333 | ierr=NF_CLOSE(outfid) |
---|
334 | if (ierr.ne.NF_NOERR) then |
---|
335 | write(*,*) 'Error, failed to close output file ',outfile |
---|
336 | endif |
---|
337 | |
---|
338 | |
---|
339 | end program |
---|