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

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

SL: small modif on angmom (Titan/Venus? Tools)

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