source: trunk/LMDZ.VENUS/Tools/angmom.F90 @ 816

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

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

File size: 25.7 KB
Line 
1program angmom
2
3! SL 12/2009:
4! This program reads 4D (lon-lat-alt-time) fields directly from LMD (or CAM) outputs
5! without regrid : histmth OR from files recast in log P coordinates (_P)
6! (+ dynzon for LMD if present, but beware of recast consistency)
7!
8! it computes:
9! dmass -- 4D -- mass of each cell
10! osam  -- 4D -- specific angular momentum (omega term)
11! rsam  -- 4D -- specific angular momentum (zonal wind term)
12! oaam  -- 1D -- integrated angular momentum (omega term)
13! raam  -- 1D -- integrated angular momentum (zonal wind term)
14! tmou  -- 1D -- Mountain torque
15! tbls  -- 1D -- Surface friction torque IF duvdf is present
16!                                       or if simple friction
17! tdyn  -- 1D -- Dynamics torque IF dudyn is present
18! tgwo  -- 1D -- Orographic GW torque IF dugwo is present
19! tgwno -- 1D -- Non-Orographic GW torque IF dugwno is present
20!
21! Minimal requirements and dependencies:
22! The dataset must include the following data:
23! - surface pressure and surface geopotential
24! - zonal wind
25! Optional: dudyn, duvdf, dugwo, dugwno (acceleration terms from physiq param)
26
27implicit none
28
29include "netcdf.inc" ! NetCDF definitions
30
31character (len=128) :: infile ! input file name
32character (len=128) :: dzfile ! input dynzon file name
33character (len=128) :: outfile ! output file name
34
35character (len=64) :: text ! to store some text
36integer infid ! NetCDF input file ID
37integer dzfid ! NetCDF input dynzon file ID
38integer outfid ! NetCDF output file ID
39integer lon_dimid,lat_dimid,alt_dimid,time_dimid ! NetCDF dimension IDs
40integer latdz_dimid,altdz_dimid,timedz_dimid ! NetCDF dimension IDs
41integer lon_varid,lat_varid,alt_varid,time_varid, tmpvarid
42integer latdz_varid,altdz_varid,timedz_varid
43integer              :: datashape1d ! shape of 1D datasets
44integer,dimension(2) :: datashape2d ! shape of 2D datasets
45integer,dimension(3) :: datashape3d ! shape of 3D datasets
46integer,dimension(4) :: datashape4d ! shape of 4D datasets
47
48real :: miss_val=9.99e+33 ! special "missing value" to specify missing data
49real,parameter :: miss_val_def=9.99e+33 ! default value for "missing value"
50real :: pi
51real,dimension(:),allocatable :: lon ! longitude
52integer lonlength ! # of grid points along longitude
53real,dimension(:),allocatable :: lat ! latitude
54real,dimension(:),allocatable :: latrad ! latitude in radian
55integer latlength ! # of grid points along latitude
56real,dimension(:),allocatable :: plev ! Pressure levels (Pa)
57integer altlength ! # of grid point along presnivs (of input datasets)
58real,dimension(:),allocatable :: time ! time
59integer timelength ! # of points along time
60real,dimension(:,:,:),allocatable :: ps ! surface pressure
61real,dimension(:,:,:),allocatable :: phis3d ! surface geopotential
62real,dimension(:,:),allocatable :: phis ! surface geopotential
63real,dimension(:,:,:,:),allocatable :: temp ! atmospheric temperature
64real,dimension(:,:,:,:),allocatable :: vitu ! zonal wind (in m/s)
65real,dimension(:,:,:,:),allocatable :: vitv ! meridional wind (in m/s)
66
67real,dimension(:,:,:,:),allocatable :: duvdf  ! Friction in BL
68real,dimension(:,:,:,:),allocatable :: dudyn  ! Dynamics
69real,dimension(:,:,:,:),allocatable :: dugwo  ! Orographic Gravity Waves
70real,dimension(:,:,:,:),allocatable :: dugwno ! Non-Orographic Gravity Waves
71
72real,dimension(:,:,:),allocatable :: dmcdyn  ! Dynamics dAM from dynzon
73real,dimension(:,:,:),allocatable :: dmcdis  ! Dissip dAM from dynzon
74real,dimension(:,:,:),allocatable :: dmcspg  ! Sponge dAM from dynzon
75real,dimension(:,:,:),allocatable :: dmcphy  ! Physics dAM from dynzon
76
77real,dimension(:,:,:,:),allocatable :: rayon ! distance to center (m)
78real,dimension(:,:,:,:),allocatable :: grav ! gravity field (m s-2)
79real,dimension(:,:,:,:),allocatable :: dmass ! mass in cell (kg)
80real,dimension(:,:,:,:),allocatable :: osam ! planetary rotation specific ang (m2/s)
81real,dimension(:,:,:,:),allocatable :: rsam ! zonal wind specific ang (m2/s)
82real,dimension(:),allocatable :: oaam ! planetary rotation total ang (kg m2/s)
83real,dimension(:),allocatable :: raam ! zonal wind total ang (kg m2/s)
84real,dimension(:),allocatable :: tmou ! mountain torque (kg m2/s2)
85real,dimension(:),allocatable :: tdyn ! dynamics torque (kg m2/s2)
86real,dimension(:),allocatable :: tbls ! friction torque (kg m2/s2)
87real,dimension(:),allocatable :: tgwo ! oro GW torque (kg m2/s2)
88real,dimension(:),allocatable :: tgwno! non-oro GW torque (kg m2/s2)
89real,dimension(:),allocatable :: tdyndz ! dynamics torque (kg m2/s2) from dynzon
90real,dimension(:),allocatable :: tdisdz ! dissip torque (kg m2/s2) from dynzon
91real,dimension(:),allocatable :: tspgdz ! sponge torque (kg m2/s2) from dynzon
92real,dimension(:),allocatable :: tphydz ! physics torque (kg m2/s2) from dynzon
93
94! for angular momentum budget normalisation
95real,parameter :: hadley=1.e18
96real,parameter :: hadday=1.e25
97
98integer ierr,ierr1,ierr2 ! NetCDF routines return codes
99integer i,j,ilon,ilat,ilev,itim ! for loops
100integer idlsurf ! for option ideal surface
101logical :: flag_duvdf,flag_dudyn,flag_dugwo,flag_dugwno,lmdflag,dzflag
102
103real :: deltalat,deltalon ! lat and lon intervals in radians
104real :: tmpp ! temporary pressure
105real :: dz ! altitude diff
106real :: signe ! orientation of lon axis for mountain torque computation
107real :: norm ! for dynzon
108
109include "planet.h"
110
111!===============================================================================
112! 1. Input parameters
113!===============================================================================
114
115pi = 2.*asin(1.)
116
117write(*,*) ""
118write(*,*) "You are working on the atmosphere of ",planet
119
120!===============================================================================
121! 1.1 Input file
122!===============================================================================
123
124write(*,*) ""
125write(*,*) " Program valid for Venus or Titan LMD, or Venus CAM output files"
126write(*,*) "Enter input file name:"
127
128read(*,'(a128)') infile
129write(*,*) ""
130
131! open input file
132
133ierr = NF_OPEN(infile,NF_NOWRITE,infid)
134if (ierr.ne.NF_NOERR) then
135   write(*,*) 'ERROR: Pb opening file ',trim(infile)
136   stop ""
137endif
138
139!===============================================================================
140! 1.2 Get grids in lon,lat,alt(pressure),time
141!===============================================================================
142
143call get_iddim(infid,lat_varid,latlength,lon_varid,lonlength,&
144                     alt_varid,altlength,time_varid,timelength,lmdflag )
145
146allocate(lon(lonlength))
147ierr=NF_GET_VAR_REAL(infid,lon_varid,lon)
148if (ierr.ne.NF_NOERR) stop "Error: Failed reading longitude"
149if(lon(1).gt.lon(2)) then
150  signe=-1.
151else
152  signe=1.
153endif
154
155allocate(lat(latlength))
156ierr=NF_GET_VAR_REAL(infid,lat_varid,lat)
157if (ierr.ne.NF_NOERR) stop "Error: Failed reading lat"
158
159allocate(latrad(latlength))
160latrad = lat*pi/180.
161
162! Lat, lon pressure intervals
163deltalat = abs(latrad(2)-latrad(1))
164deltalon = abs(lon(2)-lon(1))*pi/180.
165
166allocate(plev(altlength))
167ierr=NF_GET_VAR_REAL(infid,alt_varid,plev)
168if (ierr.ne.NF_NOERR) stop "Error: Failed reading pressure levels"
169if(.not.lmdflag) then
170! in CAM files, pressure in mbar and reversed...
171  call reverselev(altlength,plev)
172  plev=plev*100.  ! convertion to Pa
173endif
174
175allocate(time(timelength))
176ierr=NF_GET_VAR_REAL(infid,time_varid,time)
177if (ierr.ne.NF_NOERR) stop "Error: Failed reading time"
178
179! Time axis IN PLANET DAYS
180
181if(.not.lmdflag) then
182! in CAM files, time in Earth days...
183!   => seconds
184  time=time*86400.
185endif
186time=time/localday
187
188!===============================================================================
189! 1.3 dynzon file if present
190!===============================================================================
191
192! RAJOUTER UN TEST COHERENCE _P...
193
194dzflag=.false.
195
196if(lmdflag) then
197
198write(*,*) "Enter dynzon file name (<return> if not present):"
199
200read(*,'(a128)') dzfile
201write(*,*) ""
202
203if(dzfile.ne."") then
204
205! open dynzon file
206ierr = NF_OPEN(dzfile,NF_NOWRITE,dzfid)
207if (ierr.ne.NF_NOERR) then
208   write(*,*) 'ERROR: Pb opening file ',trim(dzfile)
209else
210   dzflag=.true.
211endif
212
213endif ! dzfile.ne.""
214endif ! lmdflag
215
216!===============================================================================
217! 1.4 Get output file name
218!===============================================================================
219write(*,*) ""
220!write(*,*) "Enter output file name"
221!read(*,*) outfile
222outfile=infile(1:len_trim(infile)-3)//"_GAM.nc"
223write(*,*) "Output file name is: "//trim(outfile)
224
225
226
227!===============================================================================
228! 2.1 Store needed fields
229!===============================================================================
230
231!===============================================================================
232! 2.1.1 Surface pressure and geopotential
233!===============================================================================
234allocate(ps(lonlength,latlength,timelength))
235allocate(phis3d(lonlength,latlength,timelength))
236allocate(phis(lonlength,latlength))
237
238text="PS"
239call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2)
240if (ierr1.ne.NF_NOERR) then
241  write(*,*) "  looking for psol instead... "
242  text="psol"
243  call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2)
244  if (ierr1.ne.NF_NOERR) stop "Error: Failed to get psol ID"
245endif
246if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface pressure"
247if((.not.lmdflag).and.(planet.eq."Venus")) call reverse3d(lonlength,latlength,timelength,ps)
248
249text="PHIS"
250call get_var3d(infid,lonlength,latlength,timelength,text,phis3d,ierr1,ierr2)
251if (ierr1.ne.NF_NOERR) then
252  write(*,*) " Failed to get PHIS ID (3d), looking for phis (2d) instead... "
253  text="phis"
254  call get_var2d(infid,lonlength,latlength,text,phis,ierr1,ierr2)
255  if (ierr1.ne.NF_NOERR) stop "Error: Failed to get phis ID"
256  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface geopotential"
257else
258  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface geopotential"
259  phis(:,:)=phis3d(:,:,1)
260  if((.not.lmdflag).and.(planet.eq."Venus")) call reverse2d(lonlength,latlength,phis)
261endif
262
263!===============================================================================
264! 2.1.3 Winds
265!===============================================================================
266allocate(vitu(lonlength,latlength,altlength,timelength))
267
268! zonal wind vitu / U (in m/s)
269if(lmdflag) then
270  text="vitu"
271else
272  text="U"
273endif
274
275call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitu,ierr1,ierr2)
276if (ierr1.ne.NF_NOERR) stop "Error: Failed to get U ID"
277if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zonal wind"
278
279if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,vitu)
280
281!===============================================================================
282! 2.1.4 Altitude above areoide
283!===============================================================================
284! Only needed if g(z) on Titan...
285
286!allocate(za(lonlength,latlength,altlength,timelength))
287
288!text="zareoid"
289!call get_var4d(infid,lonlength,latlength,altlength,timelength,text,za,ierr1,ierr2)
290!if (ierr1.ne.NF_NOERR) stop "Error: Failed to get za ID"
291!if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zareoid"
292
293!===============================================================================
294! 2.1.5 Friction in Boundary layer
295!===============================================================================
296allocate(duvdf(lonlength,latlength,altlength,timelength))
297
298if(lmdflag) then
299  text="duvdf"
300else
301  text="DUVDF"
302endif
303call get_var4d(infid,lonlength,latlength,altlength,timelength,text,duvdf,ierr1,ierr2)
304if (ierr1.ne.NF_NOERR) then
305  write(*,*) "Failed to get duvdf ID"
306  flag_duvdf = .false.
307
308if(.not.lmdflag) then
309! IDEAL FRICTION ?
310write(*,*) ""
311write(*,*) " Is the friction at surface ideal ? 0 for no, 1 for yes"
312write(*,*) " duvdf = -u/(86400*30) in surface layer"
313write(*,*) " timescale hard-coded: 30 Edays"
314read(*,'(i1)') idlsurf
315!write(*,*) ""
316!write(*,*) " ASSUMED YES ! "
317!idlsurf=1
318write(*,*) ""
319else
320idlsurf=0
321endif
322
323if (idlsurf.eq.1) then
324  flag_duvdf = .true.
325  duvdf = 0.
326  ilev=1
327do ilon=1,lonlength
328 do ilat=1,latlength
329  do itim=1,timelength
330   duvdf(ilon,ilat,ilev,itim) = vitu(ilon,ilat,ilev,itim) * (-1.)/(86400.*30.)
331  enddo
332 enddo
333enddo
334endif
335
336else !err1
337  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading duvdf"
338  if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,duvdf)
339  flag_duvdf = .true.
340endif !err1
341
342!===============================================================================
343! 2.1.6 Orographic and Non-orographic gravity waves
344!===============================================================================
345allocate(dugwo(lonlength,latlength,altlength,timelength))
346allocate(dugwno(lonlength,latlength,altlength,timelength))
347
348if(lmdflag) then
349
350text="dugwo"
351call get_var4d(infid,lonlength,latlength,altlength,timelength,text,dugwo,ierr1,ierr2)
352if (ierr1.ne.NF_NOERR) then
353  write(*,*) "Failed to get dugwo ID"
354  flag_dugwo = .false.
355else
356  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dugwo"
357  flag_dugwo = .true.
358endif
359
360text="dugwno"
361call get_var4d(infid,lonlength,latlength,altlength,timelength,text,dugwno,ierr1,ierr2)
362if (ierr1.ne.NF_NOERR) then
363  write(*,*) "Failed to get dugwno ID"
364  flag_dugwno = .false.
365else
366  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dugwno"
367  flag_dugwno = .true.
368endif
369
370else   ! lmdflag
371  print*,"dugwo and dugwno not in CAM simulations"
372endif  ! lmdflag
373
374!===============================================================================
375! 2.1.7 Dynamics (includes the mountain torque...)
376!===============================================================================
377allocate(dudyn(lonlength,latlength,altlength,timelength))
378
379if(lmdflag) then
380  text="dudyn"
381else
382  text="DUDYN"
383endif
384call get_var4d(infid,lonlength,latlength,altlength,timelength,text,dudyn,ierr1,ierr2)
385if (ierr1.ne.NF_NOERR) then
386  write(*,*) "Failed to get dudyn ID"
387  flag_dudyn = .false.
388else
389  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dudyn"
390  if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,dudyn)
391  flag_dudyn = .true.
392endif
393
394!===============================================================================
395! 2.1.8 Dynzon dAM/dt
396!===============================================================================
397if(dzflag) then
398
399allocate(dmcdyn(latlength-1,altlength,timelength))
400allocate(dmcdis(latlength-1,altlength,timelength))
401allocate(dmcspg(latlength-1,altlength,timelength))
402allocate(dmcphy(latlength-1,altlength,timelength))
403
404  text="dmcdyn"
405call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcdyn,ierr1,ierr2)
406if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcdyn ID"
407if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcdyn"
408
409  text="dmcdis"
410call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcdis,ierr1,ierr2)
411if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcdis ID"
412if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcdis"
413
414  text="dmcspg"
415call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcspg,ierr1,ierr2)
416if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcspg ID"
417if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcspg"
418
419  text="dmcphy"
420call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcphy,ierr1,ierr2)
421if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcphy ID"
422if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcphy"
423
424endif ! dzflag
425
426!===============================================================================
427! 2.2 Computations
428!===============================================================================
429
430!===============================================================================
431! 2.2.2 Mass in cells
432!===============================================================================
433allocate(rayon(lonlength,latlength,altlength,timelength))
434allocate(grav(lonlength,latlength,altlength,timelength))
435allocate(dmass(lonlength,latlength,altlength,timelength))
436
437do itim=1,timelength
438do ilon=1,lonlength
439 do ilat=1,latlength
440  do ilev=1,altlength
441! Need to be consistent with GCM computations
442!    if (za(ilon,ilat,ilev,itim).ne.miss_val) then
443     rayon(ilon,ilat,ilev,itim) = a0
444!     rayon(ilon,ilat,ilev,itim) = za(ilon,ilat,ilev,itim) + a0
445      grav(ilon,ilat,ilev,itim) = g0*a0*a0 &
446                 /(rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim))
447!    else
448!     rayon(ilon,ilat,ilev,itim) = miss_val
449!      grav(ilon,ilat,ilev,itim) = miss_val
450!    endif
451  enddo
452 enddo
453enddo
454enddo ! timelength
455
456call cellmass(infid,latlength,lonlength,altlength,timelength, &
457              lmdflag,deltalon,deltalat,latrad,plev,ps,grav,rayon, &
458              dmass )
459
460!===============================================================================
461! 2.2.3 Specific angular momentum
462!===============================================================================
463allocate(osam(lonlength,latlength,altlength,timelength))
464allocate(rsam(lonlength,latlength,altlength,timelength))
465
466do itim=1,timelength
467
468do ilon=1,lonlength
469 do ilat=1,latlength
470  do ilev=1,altlength
471    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
472      osam(ilon,ilat,ilev,itim) = &
473         rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim) &
474       * cos(latrad(ilat))*cos(latrad(ilat)) &
475       * omega
476      rsam(ilon,ilat,ilev,itim) = vitu(ilon,ilat,ilev,itim) &
477       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))
478    else
479      osam(ilon,ilat,ilev,itim) = miss_val
480      rsam(ilon,ilat,ilev,itim) = miss_val
481    endif
482  enddo
483 enddo
484enddo
485
486enddo ! timelength
487
488!===============================================================================
489! 2.2.4 Total angular momentum
490!===============================================================================
491allocate(oaam(timelength))
492allocate(raam(timelength))
493
494do itim=1,timelength
495
496oaam(itim) = 0.
497raam(itim) = 0.
498do ilon=1,lonlength
499 do ilat=1,latlength
500  do ilev=1,altlength
501      oaam(itim) = oaam(itim) &
502       + osam(ilon,ilat,ilev,itim)/ hadday * dmass(ilon,ilat,ilev,itim)
503      raam(itim) = raam(itim) &
504       + rsam(ilon,ilat,ilev,itim)/ hadday * dmass(ilon,ilat,ilev,itim)
505  enddo
506 enddo
507enddo
508if (oaam(itim).eq.0.) then
509  oaam(itim) = miss_val
510  raam(itim) = miss_val
511endif
512
513enddo ! timelength
514
515!===============================================================================
516! 2.2.5 Mountain, friction, dynamics and GW torques
517!===============================================================================
518allocate(tmou(timelength))
519if (flag_dudyn)  allocate(tdyn(timelength))
520if (flag_duvdf)  allocate(tbls(timelength))
521if (flag_dugwo)  allocate(tgwo(timelength))
522if (flag_dugwno) allocate(tgwno(timelength))
523
524do itim=1,timelength
525
526tmou(itim) = 0.
527do ilon=1,lonlength
528 do ilat=2,latlength-1
529      if (ilon.eq.lonlength) then
530         dz = (phis(   1,ilat)     -phis(ilon,ilat))/g0
531       tmpp = (ps(     1,ilat,itim)  +ps(ilon,ilat,itim))/2.
532      else
533         dz = (phis(ilon+1,ilat)   -phis(ilon,ilat))/g0
534       tmpp = (ps(ilon+1,ilat,itim)  +ps(ilon,ilat,itim))/2.
535      endif
536      tmou(itim) = tmou(itim) + a0*a0* deltalat*cos(latrad(ilat)) &
537       * signe*dz*tmpp &
538       / hadley
539 enddo
540enddo
541
542if (flag_dudyn) then
543tdyn(itim) = 0.
544do ilon=1,lonlength
545 do ilat=1,latlength
546  do ilev=1,altlength
547    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
548      tdyn(itim) = tdyn(itim) + dudyn(ilon,ilat,ilev,itim) &
549       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
550       * dmass(ilon,ilat,ilev,itim) &
551       / hadley
552    endif
553  enddo
554 enddo
555enddo
556  if (tdyn(itim).eq.0.) tdyn(itim) = miss_val
557endif
558
559if (flag_duvdf) then
560tbls(itim) = 0.
561do ilon=1,lonlength
562 do ilat=1,latlength
563  do ilev=1,altlength
564    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
565      tbls(itim) = tbls(itim) + duvdf(ilon,ilat,ilev,itim) &
566       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
567       * dmass(ilon,ilat,ilev,itim) &
568       / hadley
569    endif
570  enddo
571 enddo
572enddo
573  if (tbls(itim).eq.0.) tbls(itim) = miss_val
574endif
575
576if (flag_dugwo) then
577tgwo(itim) = 0.
578do ilon=1,lonlength
579 do ilat=1,latlength
580  do ilev=1,altlength
581    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
582      tgwo(itim) = tgwo(itim) + dugwo(ilon,ilat,ilev,itim) &
583       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
584       * dmass(ilon,ilat,ilev,itim) &
585       / hadley
586    endif
587  enddo
588 enddo
589enddo
590  if (tgwo(itim).eq.0.) tgwo(itim) = miss_val
591endif
592
593if (flag_dugwno) then
594tgwno(itim) = 0.
595do ilon=1,lonlength
596 do ilat=1,latlength
597  do ilev=1,altlength
598    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
599      tgwno(itim) = tgwno(itim) + dugwno(ilon,ilat,ilev,itim) &
600         * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))       &
601         * dmass(ilon,ilat,ilev,itim)  &
602       / hadley
603    endif
604  enddo
605 enddo
606enddo
607  if (tgwno(itim).eq.0.) tgwno(itim) = miss_val
608endif
609
610enddo ! timelength
611
612!===============================================================================
613! 2.2.6 Torques from dynzon
614!===============================================================================
615if(dzflag) then
616
617allocate(tdyndz(timelength))
618allocate(tdisdz(timelength))
619allocate(tspgdz(timelength))
620allocate(tphydz(timelength))
621norm=2./3.*a0*a0*omega
622
623do itim=1,timelength
624
625tdyndz(itim) = 0.
626tdisdz(itim) = 0.
627tspgdz(itim) = 0.
628tphydz(itim) = 0.
629do ilon=1,lonlength
630 do ilat=2,latlength
631  do ilev=1,altlength
632      tdyndz(itim) = tdyndz(itim) + &
633     (dmcdyn(ilat-1,ilev,itim)+dmcdyn(ilat,ilev,itim))/(2*lonlength) &
634             * norm * dmass(ilon,ilat,ilev,itim) / hadley
635      tdisdz(itim) = tdisdz(itim) + &
636     (dmcdis(ilat-1,ilev,itim)+dmcdis(ilat,ilev,itim))/(2*lonlength) &
637             * norm * dmass(ilon,ilat,ilev,itim) / hadley
638      tspgdz(itim) = tspgdz(itim) + &
639     (dmcspg(ilat-1,ilev,itim)+dmcspg(ilat,ilev,itim))/(2*lonlength) &
640             * norm * dmass(ilon,ilat,ilev,itim) / hadley
641      tphydz(itim) = tphydz(itim) + &
642     (dmcphy(ilat-1,ilev,itim)+dmcphy(ilat,ilev,itim))/(2*lonlength) &
643             * norm * dmass(ilon,ilat,ilev,itim) / hadley
644  enddo
645 enddo
646enddo
647  if (tdyndz(itim).eq.0.) tdyndz(itim) = miss_val
648  if (tdisdz(itim).eq.0.) tdisdz(itim) = miss_val
649  if (tspgdz(itim).eq.0.) tspgdz(itim) = miss_val
650  if (tphydz(itim).eq.0.) tphydz(itim) = miss_val
651
652enddo ! timelength
653
654endif ! dzflag
655
656print*,"End of computations"
657
658!===============================================================================
659! 3. Create output file
660!===============================================================================
661
662! Create output file
663ierr=NF_CREATE(outfile,NF_CLOBBER,outfid)
664if (ierr.ne.NF_NOERR) then
665  write(*,*)"Error: could not create file ",outfile
666  stop
667endif
668
669!===============================================================================
670! 3.1. Define and write dimensions
671!===============================================================================
672
673call write_dim(outfid,lonlength,latlength,altlength,timelength, &
674    lon,lat,plev,time,lon_dimid,lat_dimid,alt_dimid,time_dimid)
675
676!===============================================================================
677! 3.2. Define and write variables
678!===============================================================================
679
680! 1D Variables
681
682datashape1d=time_dimid
683 
684call write_var1d(outfid,datashape1d,timelength,&
685                "oaam      ", "Total rotation ang  ","E25kgm2s-1",miss_val,&
686                 oaam )
687
688call write_var1d(outfid,datashape1d,timelength,&
689                "raam      ", "Total wind ang      ","E25kgm2s-1",miss_val,&
690                 raam )
691
692call write_var1d(outfid,datashape1d,timelength,&
693                "tmou      ", "Mountain torque     ","E18kgm2s-2",miss_val,&
694                 tmou )
695
696if (flag_dudyn) then
697call write_var1d(outfid,datashape1d,timelength,&
698                "tdyn      ", "Dynamics torque     ","E18kgm2s-2",miss_val,&
699                 tdyn )
700endif
701
702if (flag_duvdf) then
703call write_var1d(outfid,datashape1d,timelength,&
704                "tbls      ", "Friction torque     ","E18kgm2s-2",miss_val,&
705                 tbls )
706endif
707
708if (flag_dugwo) then
709call write_var1d(outfid,datashape1d,timelength,&
710                "tgwo      ", "Orographic GW torque","E18kgm2s-2",miss_val,&
711                 tgwo )
712endif
713
714if (flag_dugwno) then
715call write_var1d(outfid,datashape1d,timelength,&
716                "tgwno     ", "Non-orogr. GW torque","E18kgm2s-2",miss_val,&
717                 tgwno )
718endif
719
720if(dzflag) then
721
722call write_var1d(outfid,datashape1d,timelength,&
723                "tdyndz    ", "Dynamics torque DZ  ","E18kgm2s-2",miss_val,&
724                 tdyndz )
725
726call write_var1d(outfid,datashape1d,timelength,&
727                "tdisdz    ", "Dissip torque DZ    ","E18kgm2s-2",miss_val,&
728                 tdisdz )
729
730call write_var1d(outfid,datashape1d,timelength,&
731                "tspgdz    ", "Sponge torque DZ    ","E18kgm2s-2",miss_val,&
732                 tspgdz )
733
734call write_var1d(outfid,datashape1d,timelength, &
735                "tphydz    ", "Physics torque DZ   ","E18kgm2s-2",miss_val,&
736                 tphydz )
737
738endif ! dzflag
739
740! 2D variables
741
742datashape2d(1)=lon_dimid
743datashape2d(2)=lat_dimid
744
745call write_var2d(outfid,datashape2d,lonlength,latlength,&
746                "phis      ", "Surface geopot      ","m2 s-2    ",miss_val,&
747                 phis )
748
749! 3D variables
750
751datashape3d(1)=lon_dimid
752datashape3d(2)=lat_dimid
753datashape3d(3)=time_dimid
754
755call write_var3d(outfid,datashape3d,lonlength,latlength,timelength,&
756                "ps        ", "Surface pressure    ","Pa        ",miss_val,&
757                 ps )
758
759! 4D variables
760
761datashape4d(1)=lon_dimid
762datashape4d(2)=lat_dimid
763datashape4d(3)=alt_dimid
764datashape4d(4)=time_dimid
765
766call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,&
767                "dmass     ", "Mass                ","kg        ",miss_val,&
768                 dmass )
769
770call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,&
771                "osam      ", "Specific rotat ang  ","E25m2s-1  ",miss_val,&
772                 osam )
773
774call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,&
775                "rsam      ", "Specific wind ang   ","E25m2s-1  ",miss_val,&
776                 rsam )
777
778!!!! Close output file
779ierr=NF_CLOSE(outfid)
780if (ierr.ne.NF_NOERR) then
781  write(*,*) 'Error, failed to close output file ',outfile
782endif
783
784
785end program
Note: See TracBrowser for help on using the repository browser.