source: trunk/LMDZ.TITAN/Tools/psi.F90 @ 819

Last change on this file since 819 was 816, checked in by slebonnois, 12 years ago

SL: tools for postprocessing (Veznus and Titan); see DOC/documentation/vt-tools.pdf

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