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

Last change on this file since 3471 was 1454, checked in by slebonnois, 9 years ago

SL: corrections in Titan and Venus tools

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