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

Last change on this file since 1379 was 1356, checked in by slebonnois, 10 years ago

SL: update to newstart/start2archive tools in Venus+Titan / additional diagnostics in radiative fluxes for Titan

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