1 | program energy |
---|
2 | |
---|
3 | ! SL 12/2009: |
---|
4 | ! This program reads 4D (lon-lat-alt-time) fields directly from LMD outputs without regrid : histmth |
---|
5 | ! |
---|
6 | ! it computes: |
---|
7 | ! dmass -- 4D -- mass of each cell |
---|
8 | ! sek -- 4D -- specific kinetic energy |
---|
9 | ! ek -- 1D -- integrated kinetic energy |
---|
10 | ! sep -- 4D -- specific potential energy |
---|
11 | ! ep -- 1D -- integrated potential energy |
---|
12 | ! |
---|
13 | ! Minimal requirements and dependencies: |
---|
14 | ! The dataset must include the following data: |
---|
15 | ! - surface pressure |
---|
16 | ! - atmospheric temperature |
---|
17 | ! - zonal and meridional winds |
---|
18 | |
---|
19 | ! VERTICAL WIND SPEED IS NEGLECTED IN KINETIC ENERGY |
---|
20 | |
---|
21 | implicit none |
---|
22 | |
---|
23 | include "netcdf.inc" ! NetCDF definitions |
---|
24 | |
---|
25 | character (len=128) :: infile ! input file name (name_P.nc) |
---|
26 | character (len=128) :: outfile ! output file name |
---|
27 | |
---|
28 | character (len=64) :: text ! to store some text |
---|
29 | integer infid ! NetCDF input file ID |
---|
30 | integer outfid ! NetCDF output file ID |
---|
31 | integer lon_dimid,lat_dimid,alt_dimid,time_dimid ! NetCDF dimension IDs |
---|
32 | integer lon_varid,lat_varid,alt_varid,time_varid |
---|
33 | integer :: datashape1d ! shape of 1D datasets |
---|
34 | integer,dimension(2) :: datashape2d ! shape of 2D datasets |
---|
35 | integer,dimension(3) :: datashape3d ! shape of 3D datasets |
---|
36 | integer,dimension(4) :: datashape4d ! shape of 4D datasets |
---|
37 | |
---|
38 | real :: miss_val=9.99e+33 ! special "missing value" to specify missing data |
---|
39 | real,parameter :: miss_val_def=9.99e+33 ! default value for "missing value" |
---|
40 | real :: pi |
---|
41 | real,dimension(:),allocatable :: lon ! longitude |
---|
42 | integer lonlength ! # of grid points along longitude |
---|
43 | real,dimension(:),allocatable :: lat ! latitude |
---|
44 | real,dimension(:),allocatable :: latrad ! latitude in radian |
---|
45 | integer latlength ! # of grid points along latitude |
---|
46 | real,dimension(:),allocatable :: plev ! Pressure levels (Pa) |
---|
47 | integer altlength ! # of grid point along presnivs (of input datasets) |
---|
48 | real,dimension(:),allocatable :: time ! time |
---|
49 | integer timelength ! # of points along time |
---|
50 | real,dimension(:,:,:),allocatable :: ps ! surface pressure |
---|
51 | real,dimension(:,:,:,:),allocatable :: temp ! atmospheric temperature |
---|
52 | real,dimension(:,:,:,:),allocatable :: vitu ! zonal wind (in m/s) |
---|
53 | real,dimension(:,:,:,:),allocatable :: vitv ! meridional wind (in m/s) |
---|
54 | |
---|
55 | real,dimension(:,:,:,:),allocatable :: rayon ! distance to center (m) |
---|
56 | real,dimension(:,:,:,:),allocatable :: grav ! gravity field (m s-2) |
---|
57 | real,dimension(:,:,:,:),allocatable :: dmass ! mass in cell (kg) |
---|
58 | real,dimension(:,:,:,:),allocatable :: sek ! specific kinetic energy |
---|
59 | real,dimension(:),allocatable :: ek ! total kinetic energy |
---|
60 | real,dimension(:,:,:,:),allocatable :: sep ! specific potential energy |
---|
61 | real,dimension(:),allocatable :: ep ! total potential energy |
---|
62 | |
---|
63 | integer ierr,ierr1,ierr2 ! NetCDF routines return codes |
---|
64 | integer i,j,ilon,ilat,ilev,itim ! for loops |
---|
65 | |
---|
66 | real :: deltalat,deltalon ! lat and lon intervals in radians |
---|
67 | real,dimension(:,:,:,:),allocatable :: deltap ! pressure thickness of each layer (Pa) |
---|
68 | real :: tmpp ! temporary pressure |
---|
69 | real :: dz ! altitude diff |
---|
70 | real :: signe ! orientation of lon axis for mountain torque computation |
---|
71 | logical :: lmdflag |
---|
72 | |
---|
73 | real :: cpdet |
---|
74 | |
---|
75 | include "planet.h" |
---|
76 | |
---|
77 | !=============================================================================== |
---|
78 | ! 1. Input parameters |
---|
79 | !=============================================================================== |
---|
80 | |
---|
81 | pi = 2.*asin(1.) |
---|
82 | |
---|
83 | write(*,*) "" |
---|
84 | write(*,*) "You are working on the atmosphere of ",planet |
---|
85 | |
---|
86 | !=============================================================================== |
---|
87 | ! 1.1 Input file |
---|
88 | !=============================================================================== |
---|
89 | |
---|
90 | write(*,*) "" |
---|
91 | write(*,*) " Program valid for Venus or Titan LMD, or Venus CAM output files" |
---|
92 | write(*,*) "Enter input file name:" |
---|
93 | |
---|
94 | read(*,'(a128)') infile |
---|
95 | write(*,*) "" |
---|
96 | |
---|
97 | ! open input file |
---|
98 | |
---|
99 | ierr = NF_OPEN(infile,NF_NOWRITE,infid) |
---|
100 | if (ierr.ne.NF_NOERR) then |
---|
101 | write(*,*) 'ERROR: Pb opening file ',trim(infile) |
---|
102 | stop "" |
---|
103 | endif |
---|
104 | |
---|
105 | !=============================================================================== |
---|
106 | ! 1.2 Get grids in lon,lat,alt(pressure),time |
---|
107 | !=============================================================================== |
---|
108 | |
---|
109 | call get_iddim(infid,lat_varid,latlength,lon_varid,lonlength,& |
---|
110 | alt_varid,altlength,time_varid,timelength,lmdflag ) |
---|
111 | |
---|
112 | allocate(lon(lonlength)) |
---|
113 | ierr=NF_GET_VAR_REAL(infid,lon_varid,lon) |
---|
114 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading longitude" |
---|
115 | if(lon(1).gt.lon(2)) then |
---|
116 | signe=-1. |
---|
117 | else |
---|
118 | signe=1. |
---|
119 | endif |
---|
120 | |
---|
121 | allocate(lat(latlength)) |
---|
122 | ierr=NF_GET_VAR_REAL(infid,lat_varid,lat) |
---|
123 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading lat" |
---|
124 | |
---|
125 | allocate(latrad(latlength)) |
---|
126 | latrad = lat*pi/180. |
---|
127 | |
---|
128 | ! Lat, lon pressure intervals |
---|
129 | deltalat = abs(latrad(2)-latrad(1)) |
---|
130 | deltalon = abs(lon(2)-lon(1))*pi/180. |
---|
131 | |
---|
132 | allocate(plev(altlength)) |
---|
133 | ierr=NF_GET_VAR_REAL(infid,alt_varid,plev) |
---|
134 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading presnivs (ie pressure levels)" |
---|
135 | if(.not.lmdflag) then |
---|
136 | ! in CAM files, pressure in mbar and reversed... |
---|
137 | call reverselev(altlength,plev) |
---|
138 | plev=plev*100. ! convertion to Pa |
---|
139 | endif |
---|
140 | |
---|
141 | allocate(time(timelength)) |
---|
142 | ierr=NF_GET_VAR_REAL(infid,time_varid,time) |
---|
143 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading time" |
---|
144 | |
---|
145 | ! Time axis IN PLANET DAYS |
---|
146 | |
---|
147 | if(.not.lmdflag) then |
---|
148 | ! in CAM files, time in Earth days... |
---|
149 | ! => seconds |
---|
150 | time=time*86400. |
---|
151 | endif |
---|
152 | time=time/localday |
---|
153 | |
---|
154 | !=============================================================================== |
---|
155 | ! 1.3 Get output file name |
---|
156 | !=============================================================================== |
---|
157 | write(*,*) "" |
---|
158 | !write(*,*) "Enter output file name" |
---|
159 | !read(*,*) outfile |
---|
160 | outfile=infile(1:len_trim(infile)-3)//"_NRG.nc" |
---|
161 | write(*,*) "Output file name is: "//trim(outfile) |
---|
162 | |
---|
163 | |
---|
164 | |
---|
165 | !=============================================================================== |
---|
166 | ! 2.1 Store needed fields |
---|
167 | !=============================================================================== |
---|
168 | |
---|
169 | !=============================================================================== |
---|
170 | ! 2.1.1 Surface pressure |
---|
171 | !=============================================================================== |
---|
172 | allocate(ps(lonlength,latlength,timelength)) |
---|
173 | |
---|
174 | text="ps" |
---|
175 | call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2) |
---|
176 | if (ierr1.ne.NF_NOERR) then |
---|
177 | write(*,*) " looking for psol instead... " |
---|
178 | text="psol" |
---|
179 | call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2) |
---|
180 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get psol ID" |
---|
181 | endif |
---|
182 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface pressure" |
---|
183 | if((.not.lmdflag).and.(planet.eq."Venus")) call reverse3d(lonlength,latlength,timelength,ps) |
---|
184 | |
---|
185 | !=============================================================================== |
---|
186 | ! 2.1.2 Atmospheric temperature |
---|
187 | !=============================================================================== |
---|
188 | allocate(temp(lonlength,latlength,altlength,timelength)) |
---|
189 | |
---|
190 | if(lmdflag) then |
---|
191 | text="temp" |
---|
192 | else |
---|
193 | text="T" |
---|
194 | endif |
---|
195 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,temp,ierr1,ierr2) |
---|
196 | if (ierr1.ne.NF_NOERR) then |
---|
197 | write(*,*) " looking for t instead... " |
---|
198 | text="t" |
---|
199 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,temp,ierr1,ierr2) |
---|
200 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get temperature ID" |
---|
201 | endif |
---|
202 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading temperature" |
---|
203 | |
---|
204 | if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,temp) |
---|
205 | |
---|
206 | !=============================================================================== |
---|
207 | ! 2.1.3 Winds |
---|
208 | !=============================================================================== |
---|
209 | allocate(vitu(lonlength,latlength,altlength,timelength)) |
---|
210 | allocate(vitv(lonlength,latlength,altlength,timelength)) |
---|
211 | |
---|
212 | ! zonal wind vitu (in m/s) |
---|
213 | if(lmdflag) then |
---|
214 | text="vitu" |
---|
215 | else |
---|
216 | text="U" |
---|
217 | endif |
---|
218 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitu,ierr1,ierr2) |
---|
219 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitu ID" |
---|
220 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zonal wind" |
---|
221 | |
---|
222 | if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,vitu) |
---|
223 | |
---|
224 | ! meridional wind vitv (in m/s) |
---|
225 | if(lmdflag) then |
---|
226 | text="vitv" |
---|
227 | else |
---|
228 | text="V" |
---|
229 | endif |
---|
230 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitv,ierr1,ierr2) |
---|
231 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitv ID" |
---|
232 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading meridional wind" |
---|
233 | |
---|
234 | if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,vitv) |
---|
235 | |
---|
236 | !=============================================================================== |
---|
237 | ! 2.1.4 Altitude above areoide |
---|
238 | !=============================================================================== |
---|
239 | ! Only needed if g(z) on Titan... |
---|
240 | |
---|
241 | !allocate(za(lonlength,latlength,altlength,timelength)) |
---|
242 | |
---|
243 | !text="zareoid" |
---|
244 | !call get_var4d(infid,lonlength,latlength,altlength,timelength,text,za,ierr1,ierr2) |
---|
245 | !if (ierr1.ne.NF_NOERR) stop "Error: Failed to get za ID" |
---|
246 | !if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zareoid" |
---|
247 | |
---|
248 | !=============================================================================== |
---|
249 | ! 2.2 Computations |
---|
250 | !=============================================================================== |
---|
251 | |
---|
252 | !=============================================================================== |
---|
253 | ! 2.2.2 Mass in cells |
---|
254 | !=============================================================================== |
---|
255 | allocate(rayon(lonlength,latlength,altlength,timelength)) |
---|
256 | allocate(grav(lonlength,latlength,altlength,timelength)) |
---|
257 | allocate(dmass(lonlength,latlength,altlength,timelength)) |
---|
258 | |
---|
259 | do itim=1,timelength |
---|
260 | do ilon=1,lonlength |
---|
261 | do ilat=1,latlength |
---|
262 | do ilev=1,altlength |
---|
263 | ! Need to be consistent with GCM computations |
---|
264 | ! if (za(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
265 | rayon(ilon,ilat,ilev,itim) = a0 |
---|
266 | ! rayon(ilon,ilat,ilev,itim) = za(ilon,ilat,ilev,itim) + a0 |
---|
267 | grav(ilon,ilat,ilev,itim) = g0*a0*a0 & |
---|
268 | /(rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim)) |
---|
269 | ! else |
---|
270 | ! rayon(ilon,ilat,ilev,itim) = miss_val |
---|
271 | ! grav(ilon,ilat,ilev,itim) = miss_val |
---|
272 | ! endif |
---|
273 | enddo |
---|
274 | enddo |
---|
275 | enddo |
---|
276 | enddo ! timelength |
---|
277 | |
---|
278 | call cellmass(infid,latlength,lonlength,altlength,timelength, & |
---|
279 | lmdflag,deltalon,deltalat,latrad,plev,ps,grav,rayon, & |
---|
280 | dmass ) |
---|
281 | |
---|
282 | !=============================================================================== |
---|
283 | ! 2.2.6 Specific energies |
---|
284 | !=============================================================================== |
---|
285 | allocate(sek(lonlength,latlength,altlength,timelength)) |
---|
286 | allocate(sep(lonlength,latlength,altlength,timelength)) |
---|
287 | |
---|
288 | do itim=1,timelength |
---|
289 | |
---|
290 | do ilon=1,lonlength |
---|
291 | do ilat=1,latlength |
---|
292 | do ilev=1,altlength |
---|
293 | if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
294 | if ((vitu(ilon,ilat,ilev,itim).lt.miss_val) & |
---|
295 | .and.(vitv(ilon,ilat,ilev,itim).lt.miss_val)) then |
---|
296 | sek(ilon,ilat,ilev,itim) = 0.5 * & |
---|
297 | ( vitu(ilon,ilat,ilev,itim)*vitu(ilon,ilat,ilev,itim) & |
---|
298 | + vitv(ilon,ilat,ilev,itim)*vitv(ilon,ilat,ilev,itim) ) |
---|
299 | else |
---|
300 | sek(ilon,ilat,ilev,itim) = miss_val |
---|
301 | endif |
---|
302 | if (temp(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
303 | sep(ilon,ilat,ilev,itim) = temp(ilon,ilat,ilev,itim) & |
---|
304 | * cpdet(temp(ilon,ilat,ilev,itim)) |
---|
305 | else |
---|
306 | sep(ilon,ilat,ilev,itim) = miss_val |
---|
307 | endif |
---|
308 | else |
---|
309 | sek(ilon,ilat,ilev,itim) = miss_val |
---|
310 | sep(ilon,ilat,ilev,itim) = miss_val |
---|
311 | endif |
---|
312 | enddo |
---|
313 | enddo |
---|
314 | enddo |
---|
315 | |
---|
316 | enddo ! timelength |
---|
317 | |
---|
318 | !=============================================================================== |
---|
319 | ! 2.2.7 Total energies |
---|
320 | !=============================================================================== |
---|
321 | allocate(ek(timelength)) |
---|
322 | allocate(ep(timelength)) |
---|
323 | |
---|
324 | do itim=1,timelength |
---|
325 | |
---|
326 | ek(itim) = 0. |
---|
327 | ep(itim) = 0. |
---|
328 | do ilon=1,lonlength |
---|
329 | do ilat=1,latlength |
---|
330 | do ilev=1,altlength |
---|
331 | if (sek(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
332 | ek(itim) = ek(itim) & |
---|
333 | + sek(ilon,ilat,ilev,itim) * dmass(ilon,ilat,ilev,itim) |
---|
334 | endif |
---|
335 | if (sep(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
336 | ep(itim) = ep(itim) & |
---|
337 | + sep(ilon,ilat,ilev,itim) * dmass(ilon,ilat,ilev,itim) |
---|
338 | endif |
---|
339 | enddo |
---|
340 | enddo |
---|
341 | enddo |
---|
342 | if (ek(itim).eq.0.) then |
---|
343 | ek(itim) = miss_val |
---|
344 | ep(itim) = miss_val |
---|
345 | endif |
---|
346 | |
---|
347 | enddo ! timelength |
---|
348 | |
---|
349 | print*,"End of computations" |
---|
350 | |
---|
351 | !=============================================================================== |
---|
352 | ! 3. Create output file |
---|
353 | !=============================================================================== |
---|
354 | |
---|
355 | ! Create output file |
---|
356 | ierr=NF_CREATE(outfile,NF_CLOBBER,outfid) |
---|
357 | if (ierr.ne.NF_NOERR) then |
---|
358 | write(*,*)"Error: could not create file ",outfile |
---|
359 | stop |
---|
360 | endif |
---|
361 | |
---|
362 | !=============================================================================== |
---|
363 | ! 3.1. Define and write dimensions |
---|
364 | !=============================================================================== |
---|
365 | |
---|
366 | call write_dim(outfid,lonlength,latlength,altlength,timelength, & |
---|
367 | lon,lat,plev,time,lon_dimid,lat_dimid,alt_dimid,time_dimid) |
---|
368 | |
---|
369 | !=============================================================================== |
---|
370 | ! 3.2. Define and write variables |
---|
371 | !=============================================================================== |
---|
372 | |
---|
373 | ! 1D Variables |
---|
374 | |
---|
375 | datashape1d=time_dimid |
---|
376 | |
---|
377 | call write_var1d(outfid,datashape1d,timelength,& |
---|
378 | "ek ", "Total kinetic energy","J ",miss_val,& |
---|
379 | ek ) |
---|
380 | |
---|
381 | call write_var1d(outfid,datashape1d,timelength,& |
---|
382 | "ep ", "Total pot energy ","J ",miss_val,& |
---|
383 | ep ) |
---|
384 | |
---|
385 | ! 3D variables |
---|
386 | |
---|
387 | datashape3d(1)=lon_dimid |
---|
388 | datashape3d(2)=lat_dimid |
---|
389 | datashape3d(3)=time_dimid |
---|
390 | |
---|
391 | call write_var3d(outfid,datashape3d,lonlength,latlength,timelength,& |
---|
392 | "ps ", "Surface pressure ","Pa ",miss_val,& |
---|
393 | ps ) |
---|
394 | |
---|
395 | ! 4D variables |
---|
396 | |
---|
397 | datashape4d(1)=lon_dimid |
---|
398 | datashape4d(2)=lat_dimid |
---|
399 | datashape4d(3)=alt_dimid |
---|
400 | datashape4d(4)=time_dimid |
---|
401 | |
---|
402 | call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
403 | "dmass ", "Mass ","kg ",miss_val,& |
---|
404 | dmass ) |
---|
405 | |
---|
406 | call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
407 | "sek ", "Specific kin energy ","J/kg ",miss_val,& |
---|
408 | sek ) |
---|
409 | |
---|
410 | call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
411 | "sep ", "Specific pot energy ","J/kg ",miss_val,& |
---|
412 | sep ) |
---|
413 | |
---|
414 | |
---|
415 | !!!! Close output file |
---|
416 | ierr=NF_CLOSE(outfid) |
---|
417 | if (ierr.ne.NF_NOERR) then |
---|
418 | write(*,*) 'Error, failed to close output file ',outfile |
---|
419 | endif |
---|
420 | |
---|
421 | |
---|
422 | end program |
---|