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

Last change on this file since 819 was 819, checked in by slebonnois, 12 years ago
  • newstart-VT.F and start2archive-VT.F added in dyn3d/, their name indicating that they are tools for Venus and Titan.
  • doc for these tools added
  • correction of a bug on angmom (Venus and Titan postprocessing tools)
  • one tool (lat_torques) removed because not updated until further needed
  • correction of a bug in phyetat0 (Venus)
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
238if(lmdflag) then
239text="psol"
240else
241text="PS"
242endif
243call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2)
244if (ierr1.ne.NF_NOERR) then
245  write(*,*) "  looking for ps instead... "
246  text="ps"
247  call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2)
248  if (ierr1.ne.NF_NOERR) stop "Error: Failed to get ps ID"
249endif
250if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface pressure"
251if((.not.lmdflag).and.(planet.eq."Venus")) call reverse3d(lonlength,latlength,timelength,ps)
252
253text="PHIS"
254call get_var3d(infid,lonlength,latlength,timelength,text,phis3d,ierr1,ierr2)
255if (ierr1.ne.NF_NOERR) then
256  write(*,*) " Failed to get PHIS ID (3d), looking for phis (2d) instead... "
257  text="phis"
258  call get_var2d(infid,lonlength,latlength,text,phis,ierr1,ierr2)
259  if (ierr1.ne.NF_NOERR) stop "Error: Failed to get phis ID"
260  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface geopotential"
261else
262  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface geopotential"
263  phis(:,:)=phis3d(:,:,1)
264  if((.not.lmdflag).and.(planet.eq."Venus")) call reverse2d(lonlength,latlength,phis)
265endif
266
267!===============================================================================
268! 2.1.3 Winds
269!===============================================================================
270allocate(vitu(lonlength,latlength,altlength,timelength))
271
272! zonal wind vitu / U (in m/s)
273if(lmdflag) then
274  text="vitu"
275else
276  text="U"
277endif
278
279call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitu,ierr1,ierr2)
280if (ierr1.ne.NF_NOERR) stop "Error: Failed to get U ID"
281if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zonal wind"
282
283if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,vitu)
284
285!===============================================================================
286! 2.1.4 Altitude above areoide
287!===============================================================================
288! Only needed if g(z) on Titan...
289
290!allocate(za(lonlength,latlength,altlength,timelength))
291
292!text="zareoid"
293!call get_var4d(infid,lonlength,latlength,altlength,timelength,text,za,ierr1,ierr2)
294!if (ierr1.ne.NF_NOERR) stop "Error: Failed to get za ID"
295!if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zareoid"
296
297!===============================================================================
298! 2.1.5 Friction in Boundary layer
299!===============================================================================
300allocate(duvdf(lonlength,latlength,altlength,timelength))
301
302if(lmdflag) then
303  text="duvdf"
304else
305  text="DUVDF"
306endif
307call get_var4d(infid,lonlength,latlength,altlength,timelength,text,duvdf,ierr1,ierr2)
308if (ierr1.ne.NF_NOERR) then
309  write(*,*) "Failed to get duvdf ID"
310  flag_duvdf = .false.
311
312if(.not.lmdflag) then
313! IDEAL FRICTION ?
314write(*,*) ""
315write(*,*) " Is the friction at surface ideal ? 0 for no, 1 for yes"
316write(*,*) " duvdf = -u/(86400*30) in surface layer"
317write(*,*) " timescale hard-coded: 30 Edays"
318read(*,'(i1)') idlsurf
319!write(*,*) ""
320!write(*,*) " ASSUMED YES ! "
321!idlsurf=1
322write(*,*) ""
323else
324idlsurf=0
325endif
326
327if (idlsurf.eq.1) then
328  flag_duvdf = .true.
329  duvdf = 0.
330  ilev=1
331do ilon=1,lonlength
332 do ilat=1,latlength
333  do itim=1,timelength
334   duvdf(ilon,ilat,ilev,itim) = vitu(ilon,ilat,ilev,itim) * (-1.)/(86400.*30.)
335  enddo
336 enddo
337enddo
338endif
339
340else !err1
341  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading duvdf"
342  if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,duvdf)
343  flag_duvdf = .true.
344endif !err1
345
346!===============================================================================
347! 2.1.6 Orographic and Non-orographic gravity waves
348!===============================================================================
349allocate(dugwo(lonlength,latlength,altlength,timelength))
350allocate(dugwno(lonlength,latlength,altlength,timelength))
351
352if(lmdflag) then
353
354text="dugwo"
355call get_var4d(infid,lonlength,latlength,altlength,timelength,text,dugwo,ierr1,ierr2)
356if (ierr1.ne.NF_NOERR) then
357  write(*,*) "Failed to get dugwo ID"
358  flag_dugwo = .false.
359else
360  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dugwo"
361  flag_dugwo = .true.
362endif
363
364text="dugwno"
365call get_var4d(infid,lonlength,latlength,altlength,timelength,text,dugwno,ierr1,ierr2)
366if (ierr1.ne.NF_NOERR) then
367  write(*,*) "Failed to get dugwno ID"
368  flag_dugwno = .false.
369else
370  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dugwno"
371  flag_dugwno = .true.
372endif
373
374else   ! lmdflag
375  print*,"dugwo and dugwno not in CAM simulations"
376endif  ! lmdflag
377
378!===============================================================================
379! 2.1.7 Dynamics (includes the mountain torque...)
380!===============================================================================
381allocate(dudyn(lonlength,latlength,altlength,timelength))
382
383if(lmdflag) then
384  text="dudyn"
385else
386  text="DUDYN"
387endif
388call get_var4d(infid,lonlength,latlength,altlength,timelength,text,dudyn,ierr1,ierr2)
389if (ierr1.ne.NF_NOERR) then
390  write(*,*) "Failed to get dudyn ID"
391  flag_dudyn = .false.
392else
393  if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dudyn"
394  if(.not.lmdflag) call reverse4d(lonlength,latlength,altlength,timelength,dudyn)
395  flag_dudyn = .true.
396endif
397
398!===============================================================================
399! 2.1.8 Dynzon dAM/dt
400!===============================================================================
401if(dzflag) then
402
403allocate(dmcdyn(latlength-1,altlength,timelength))
404allocate(dmcdis(latlength-1,altlength,timelength))
405allocate(dmcspg(latlength-1,altlength,timelength))
406allocate(dmcphy(latlength-1,altlength,timelength))
407
408  text="dmcdyn"
409call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcdyn,ierr1,ierr2)
410if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcdyn ID"
411if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcdyn"
412
413  text="dmcdis"
414call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcdis,ierr1,ierr2)
415if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcdis ID"
416if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcdis"
417
418  text="dmcspg"
419call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcspg,ierr1,ierr2)
420if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcspg ID"
421if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcspg"
422
423  text="dmcphy"
424call get_var3d(dzfid,latlength-1,altlength,timelength,text,dmcphy,ierr1,ierr2)
425if (ierr1.ne.NF_NOERR) stop "Error: Failed to get dmcphy ID"
426if (ierr2.ne.NF_NOERR) stop "Error: Failed reading dmcphy"
427
428endif ! dzflag
429
430!===============================================================================
431! 2.2 Computations
432!===============================================================================
433
434!===============================================================================
435! 2.2.2 Mass in cells
436!===============================================================================
437allocate(rayon(lonlength,latlength,altlength,timelength))
438allocate(grav(lonlength,latlength,altlength,timelength))
439allocate(dmass(lonlength,latlength,altlength,timelength))
440
441do itim=1,timelength
442do ilon=1,lonlength
443 do ilat=1,latlength
444  do ilev=1,altlength
445! Need to be consistent with GCM computations
446!    if (za(ilon,ilat,ilev,itim).ne.miss_val) then
447     rayon(ilon,ilat,ilev,itim) = a0
448!     rayon(ilon,ilat,ilev,itim) = za(ilon,ilat,ilev,itim) + a0
449      grav(ilon,ilat,ilev,itim) = g0*a0*a0 &
450                 /(rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim))
451!    else
452!     rayon(ilon,ilat,ilev,itim) = miss_val
453!      grav(ilon,ilat,ilev,itim) = miss_val
454!    endif
455  enddo
456 enddo
457enddo
458enddo ! timelength
459
460call cellmass(infid,latlength,lonlength,altlength,timelength, &
461              lmdflag,deltalon,deltalat,latrad,plev,ps,grav,rayon, &
462              dmass )
463
464!===============================================================================
465! 2.2.3 Specific angular momentum
466!===============================================================================
467allocate(osam(lonlength,latlength,altlength,timelength))
468allocate(rsam(lonlength,latlength,altlength,timelength))
469
470do itim=1,timelength
471
472do ilon=1,lonlength
473 do ilat=1,latlength
474  do ilev=1,altlength
475    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
476      osam(ilon,ilat,ilev,itim) = &
477         rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim) &
478       * cos(latrad(ilat))*cos(latrad(ilat)) &
479       * omega
480      rsam(ilon,ilat,ilev,itim) = vitu(ilon,ilat,ilev,itim) &
481       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))
482    else
483      osam(ilon,ilat,ilev,itim) = miss_val
484      rsam(ilon,ilat,ilev,itim) = miss_val
485    endif
486  enddo
487 enddo
488enddo
489
490enddo ! timelength
491
492!===============================================================================
493! 2.2.4 Total angular momentum
494!===============================================================================
495allocate(oaam(timelength))
496allocate(raam(timelength))
497
498do itim=1,timelength
499
500oaam(itim) = 0.
501raam(itim) = 0.
502do ilon=1,lonlength
503 do ilat=1,latlength
504  do ilev=1,altlength
505      oaam(itim) = oaam(itim) &
506       + osam(ilon,ilat,ilev,itim)/ hadday * dmass(ilon,ilat,ilev,itim)
507      raam(itim) = raam(itim) &
508       + rsam(ilon,ilat,ilev,itim)/ hadday * dmass(ilon,ilat,ilev,itim)
509  enddo
510 enddo
511enddo
512if (oaam(itim).eq.0.) then
513  oaam(itim) = miss_val
514  raam(itim) = miss_val
515endif
516
517enddo ! timelength
518
519!===============================================================================
520! 2.2.5 Mountain, friction, dynamics and GW torques
521!===============================================================================
522allocate(tmou(timelength))
523if (flag_dudyn)  allocate(tdyn(timelength))
524if (flag_duvdf)  allocate(tbls(timelength))
525if (flag_dugwo)  allocate(tgwo(timelength))
526if (flag_dugwno) allocate(tgwno(timelength))
527
528do itim=1,timelength
529
530tmou(itim) = 0.
531do ilon=1,lonlength
532 do ilat=2,latlength-1
533      if (ilon.eq.lonlength) then
534         dz = (phis(   1,ilat)     -phis(ilon,ilat))/g0
535       tmpp = (ps(     1,ilat,itim)  +ps(ilon,ilat,itim))/2.
536      else
537         dz = (phis(ilon+1,ilat)   -phis(ilon,ilat))/g0
538       tmpp = (ps(ilon+1,ilat,itim)  +ps(ilon,ilat,itim))/2.
539      endif
540      tmou(itim) = tmou(itim) + a0*a0* deltalat*cos(latrad(ilat)) &
541       * signe*dz*tmpp &
542       / hadley
543 enddo
544enddo
545
546if (flag_dudyn) then
547tdyn(itim) = 0.
548do ilon=1,lonlength
549 do ilat=1,latlength
550  do ilev=1,altlength
551    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
552      tdyn(itim) = tdyn(itim) + dudyn(ilon,ilat,ilev,itim) &
553       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
554       * dmass(ilon,ilat,ilev,itim) &
555       / hadley
556    endif
557  enddo
558 enddo
559enddo
560  if (tdyn(itim).eq.0.) tdyn(itim) = miss_val
561endif
562
563if (flag_duvdf) then
564tbls(itim) = 0.
565do ilon=1,lonlength
566 do ilat=1,latlength
567  do ilev=1,altlength
568    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
569      tbls(itim) = tbls(itim) + duvdf(ilon,ilat,ilev,itim) &
570       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
571       * dmass(ilon,ilat,ilev,itim) &
572       / hadley
573    endif
574  enddo
575 enddo
576enddo
577  if (tbls(itim).eq.0.) tbls(itim) = miss_val
578endif
579
580if (flag_dugwo) then
581tgwo(itim) = 0.
582do ilon=1,lonlength
583 do ilat=1,latlength
584  do ilev=1,altlength
585    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
586      tgwo(itim) = tgwo(itim) + dugwo(ilon,ilat,ilev,itim) &
587       * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))      &
588       * dmass(ilon,ilat,ilev,itim) &
589       / hadley
590    endif
591  enddo
592 enddo
593enddo
594  if (tgwo(itim).eq.0.) tgwo(itim) = miss_val
595endif
596
597if (flag_dugwno) then
598tgwno(itim) = 0.
599do ilon=1,lonlength
600 do ilat=1,latlength
601  do ilev=1,altlength
602    if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then
603      tgwno(itim) = tgwno(itim) + dugwno(ilon,ilat,ilev,itim) &
604         * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))       &
605         * dmass(ilon,ilat,ilev,itim)  &
606       / hadley
607    endif
608  enddo
609 enddo
610enddo
611  if (tgwno(itim).eq.0.) tgwno(itim) = miss_val
612endif
613
614enddo ! timelength
615
616!===============================================================================
617! 2.2.6 Torques from dynzon
618!===============================================================================
619if(dzflag) then
620
621allocate(tdyndz(timelength))
622allocate(tdisdz(timelength))
623allocate(tspgdz(timelength))
624allocate(tphydz(timelength))
625norm=2./3.*a0*a0*omega
626
627do itim=1,timelength
628
629tdyndz(itim) = 0.
630tdisdz(itim) = 0.
631tspgdz(itim) = 0.
632tphydz(itim) = 0.
633do ilon=1,lonlength
634 do ilat=2,latlength
635  do ilev=1,altlength
636      tdyndz(itim) = tdyndz(itim) + &
637     (dmcdyn(ilat-1,ilev,itim)+dmcdyn(ilat,ilev,itim))/(2*lonlength) &
638             * norm * dmass(ilon,ilat,ilev,itim) / hadley
639      tdisdz(itim) = tdisdz(itim) + &
640     (dmcdis(ilat-1,ilev,itim)+dmcdis(ilat,ilev,itim))/(2*lonlength) &
641             * norm * dmass(ilon,ilat,ilev,itim) / hadley
642      tspgdz(itim) = tspgdz(itim) + &
643     (dmcspg(ilat-1,ilev,itim)+dmcspg(ilat,ilev,itim))/(2*lonlength) &
644             * norm * dmass(ilon,ilat,ilev,itim) / hadley
645      tphydz(itim) = tphydz(itim) + &
646     (dmcphy(ilat-1,ilev,itim)+dmcphy(ilat,ilev,itim))/(2*lonlength) &
647             * norm * dmass(ilon,ilat,ilev,itim) / hadley
648  enddo
649 enddo
650enddo
651  if (tdyndz(itim).eq.0.) tdyndz(itim) = miss_val
652  if (tdisdz(itim).eq.0.) tdisdz(itim) = miss_val
653  if (tspgdz(itim).eq.0.) tspgdz(itim) = miss_val
654  if (tphydz(itim).eq.0.) tphydz(itim) = miss_val
655
656enddo ! timelength
657
658endif ! dzflag
659
660print*,"End of computations"
661
662!===============================================================================
663! 3. Create output file
664!===============================================================================
665
666! Create output file
667ierr=NF_CREATE(outfile,NF_CLOBBER,outfid)
668if (ierr.ne.NF_NOERR) then
669  write(*,*)"Error: could not create file ",outfile
670  stop
671endif
672
673!===============================================================================
674! 3.1. Define and write dimensions
675!===============================================================================
676
677call write_dim(outfid,lonlength,latlength,altlength,timelength, &
678    lon,lat,plev,time,lon_dimid,lat_dimid,alt_dimid,time_dimid)
679
680!===============================================================================
681! 3.2. Define and write variables
682!===============================================================================
683
684! 1D Variables
685
686datashape1d=time_dimid
687 
688call write_var1d(outfid,datashape1d,timelength,&
689                "oaam      ", "Total rotation ang  ","E25kgm2s-1",miss_val,&
690                 oaam )
691
692call write_var1d(outfid,datashape1d,timelength,&
693                "raam      ", "Total wind ang      ","E25kgm2s-1",miss_val,&
694                 raam )
695
696call write_var1d(outfid,datashape1d,timelength,&
697                "tmou      ", "Mountain torque     ","E18kgm2s-2",miss_val,&
698                 tmou )
699
700if (flag_dudyn) then
701call write_var1d(outfid,datashape1d,timelength,&
702                "tdyn      ", "Dynamics torque     ","E18kgm2s-2",miss_val,&
703                 tdyn )
704endif
705
706if (flag_duvdf) then
707call write_var1d(outfid,datashape1d,timelength,&
708                "tbls      ", "Friction torque     ","E18kgm2s-2",miss_val,&
709                 tbls )
710endif
711
712if (flag_dugwo) then
713call write_var1d(outfid,datashape1d,timelength,&
714                "tgwo      ", "Orographic GW torque","E18kgm2s-2",miss_val,&
715                 tgwo )
716endif
717
718if (flag_dugwno) then
719call write_var1d(outfid,datashape1d,timelength,&
720                "tgwno     ", "Non-orogr. GW torque","E18kgm2s-2",miss_val,&
721                 tgwno )
722endif
723
724if(dzflag) then
725
726call write_var1d(outfid,datashape1d,timelength,&
727                "tdyndz    ", "Dynamics torque DZ  ","E18kgm2s-2",miss_val,&
728                 tdyndz )
729
730call write_var1d(outfid,datashape1d,timelength,&
731                "tdisdz    ", "Dissip torque DZ    ","E18kgm2s-2",miss_val,&
732                 tdisdz )
733
734call write_var1d(outfid,datashape1d,timelength,&
735                "tspgdz    ", "Sponge torque DZ    ","E18kgm2s-2",miss_val,&
736                 tspgdz )
737
738call write_var1d(outfid,datashape1d,timelength, &
739                "tphydz    ", "Physics torque DZ   ","E18kgm2s-2",miss_val,&
740                 tphydz )
741
742endif ! dzflag
743
744! 2D variables
745
746datashape2d(1)=lon_dimid
747datashape2d(2)=lat_dimid
748
749call write_var2d(outfid,datashape2d,lonlength,latlength,&
750                "phis      ", "Surface geopot      ","m2 s-2    ",miss_val,&
751                 phis )
752
753! 3D variables
754
755datashape3d(1)=lon_dimid
756datashape3d(2)=lat_dimid
757datashape3d(3)=time_dimid
758
759call write_var3d(outfid,datashape3d,lonlength,latlength,timelength,&
760                "ps        ", "Surface pressure    ","Pa        ",miss_val,&
761                 ps )
762
763! 4D variables
764
765datashape4d(1)=lon_dimid
766datashape4d(2)=lat_dimid
767datashape4d(3)=alt_dimid
768datashape4d(4)=time_dimid
769
770call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,&
771                "dmass     ", "Mass                ","kg        ",miss_val,&
772                 dmass )
773
774call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,&
775                "osam      ", "Specific rotat ang  ","E25m2s-1  ",miss_val,&
776                 osam )
777
778call write_var4d(outfid,datashape4d,lonlength,latlength,altlength,timelength,&
779                "rsam      ", "Specific wind ang   ","E25m2s-1  ",miss_val,&
780                 rsam )
781
782!!!! Close output file
783ierr=NF_CLOSE(outfid)
784if (ierr.ne.NF_NOERR) then
785  write(*,*) 'Error, failed to close output file ',outfile
786endif
787
788
789end program
Note: See TracBrowser for help on using the repository browser.