source: trunk/UTIL/SPECTRA/spectra_analysis.f90 @ 3026

Last change on this file since 3026 was 1426, checked in by milmd, 10 years ago

UTIL/SPECTRA. Correction of some bugs with spectra average. Now tested on Ciclad.

File size: 19.0 KB
Line 
1!M. Indurain & E. Millour
2!March 2015
3
4PROGRAM spectra_analysis
5
6!Compute energy spectrum of a 2d wind field with spherepack routines.
7!v = sum(n=1..nlat-1) 0.5*b(0,n)*B(0,n)+0.5*c(0,n)*C(0,n) + sum(m=1..mmax-1)sum(n=m..nlat-1) b(m,n)*B(m,n)+c(m,n)*C(m,n)
8!where B(m,n) and C(m,n) are vector spherical harmonics
9!Then energy spectra is defined as:
10!E(n) = 0.5*|b(0,n)|^2 + 0.5*|c(0,n)|^2 + sum(m=1..n) |b(m,n)|^2 + |c(m,n)|^2
11
12!Wind field read from netcdf file follows lmdz convention:
13!  u(nlong+1,nlat) zonal eastward wind, v(nlong+1,nlat) latitudinal wind
14!Wind field needed for spherepack routines must be like:
15!  u(nlat,nlong) zonal eastward wind, v(nlong,nlat) colatitudinal wind (from North pole)
16
17use netcdf
18
19implicit none
20
21!Physical grid
22integer :: nlat,nlong,nlongp1,nalt,ntime
23double precision, dimension(:,:), allocatable :: merid_wind,zonal_wind !spherepack convention
24double precision, dimension(:,:), allocatable :: merid_wind_lmdz,zonal_wind_lmdz !lmdz convention
25!Spectral grid
26integer :: mmax !spectra truncation mmax=min(nalt,nlong/2) or min(nalt,(nlong+1)/2)
27double precision, dimension (:,:), allocatable :: br,bi,cr,ci !spectra coefficients
28double precision, dimension (:,:), allocatable :: norm,norm_b,norm_c !norm of spectra coefficients
29double precision, dimension (:), allocatable :: spectra,spectra_b,spectra_c
30double precision, dimension (:,:), allocatable :: spectra_config,spectra_b_config,spectra_c_config
31double precision, dimension (:,:), allocatable :: spectra_global,spectra_b_global,spectra_c_global
32!Spherepack declarations
33integer :: ierror,lvhaec,ldwork,lwork,l1,l2,nt=1,ityp=0,idvw,jdvw,mdab,ndab
34double precision, dimension (:), allocatable :: wvhaec,work
35double precision, dimension (:), allocatable :: dwork
36!Netcdf stuff
37integer :: idfile,idmerid_wind,idzonal_wind
38integer :: idlat,idlong,idalt,idtime
39!Input arguments
40character (len=100), dimension(:), allocatable :: arg
41character (len=100) :: file_netcdf,file_spectra
42integer :: n_indt,n_indz,meant,meanz
43integer, dimension(:), allocatable :: tab_indt,tab_indz
44character (len=100) :: name_u,name_v,name_lat,name_long,name_alt,name_time
45character (len=100) :: tab_name_u(4)=(/'u','U','vitu','zonal_wind'/)
46character (len=100) :: tab_name_v(4)=(/'v','V','vitv','merid_wind'/)
47character (len=100) :: tab_name_lat(2)=(/'latitude','lat'/)
48character (len=100) :: tab_name_long(2)=(/'longitude','lon'/)
49character (len=100) :: tab_name_alt(3)=(/'altitude','alt','presnivs'/)
50character (len=100) :: tab_name_time(3)=(/'time','Time','time_counter'/)
51logical :: is_div_rot
52integer :: narg
53!Other
54integer :: i,j,t,z,mt,mz
55integer :: number_config        !number of (t,z) configurations wanted by user (from 1 to n_indt*n_indz))
56integer :: number_local         !number of (t,z) sub-configurations to mean (from 1 to (meant+1)*(meanz+1))
57integer :: number_global        !number of total (t,z) configurations to compute (from 1 to n_indt*n_indz*(meant+1)*(meanz+1))
58logical :: is_file,first_reading=.true.
59integer, dimension(:), allocatable :: tmp_tab_int
60character (len=100) :: tmp_char
61
62!**********
63!Initialisation
64file_netcdf = '??@@??'
65file_spectra = 'spectra'
66n_indt = 0
67n_indz = 0
68allocate(tab_indt(20))
69allocate(tab_indz(20))
70allocate(tmp_tab_int(20))
71tab_indt(:) = 1
72tab_indz(:) = 1
73meant = 0
74meanz = 0
75name_u = '??@@??'
76name_v = '??@@??'
77name_lat = '??@@??'
78name_long = '??@@??'
79name_alt = '??@@??'
80name_time = '??@@??'
81is_div_rot=.false.
82
83!**********
84!Input reading
85narg = command_argument_count()
86allocate(arg(narg))
87do i = 1, narg
88  call get_command_argument(i,arg(i))
89end do
90i = 1
91do while (i .le. narg)
92  if (arg(i) == '-o' .or. arg(i) == '-t' .or. arg(i) == '-z' .or. arg(i) == '-mt' .or. arg(i) == '-mz'&
93                .or. arg(i) == '-u' .or. arg(i) == '-v' .or. arg(i) == '-lat' .or. arg(i) == '-lon' &
94                .or. arg(i) == '-alt' .or. arg(i) == '-time') then
95    select case(arg(i))
96    case('-t')
97      n_indt = n_indt + 1
98      read(arg(i+1),'(i10)' ) tab_indt(n_indt) !temporal indice                 
99    case('-z')
100      n_indz = n_indz + 1
101      read(arg(i+1),'(i10)' ) tab_indz(n_indz) !vertical indice
102    case('-mt')
103      read(arg(i+1),'(i10)' ) meant !extent of temporal mean
104    case('-mz')
105      read(arg(i+1),'(i10)' ) meanz !extent of vertical mean
106    case('-o')
107      file_spectra=arg(i+1) !spectra file
108    case('-u')
109      read(arg(i+1),'(a100)' ) name_u !name of zonal wind
110    case('-v')
111      read(arg(i+1),'(a100)' ) name_v !name of meridional wind
112    case('-lat')
113      read(arg(i+1),'(a100)' ) name_lat !name of latitude axis
114    case('-lon')
115      read(arg(i+1),'(a100)' ) name_long !name of longitude axis
116    case('-alt')
117      read(arg(i+1),'(a100)' ) name_alt !name of altitude axis
118    case('-time')
119      read(arg(i+1),'(a100)' ) name_time !name of time axis
120    end select
121    i = i + 2
122    elseif (arg(i) == '-divrot') then
123      is_div_rot = .true.
124      i = i + 1
125    elseif (arg(i) == '-h' .or. arg(i) == '--help') then
126      print*,'Usage'
127      print*,'spectra_analysis netcdfFile [option]'
128      print*,'[-h or --help]    : brief help'
129      print*,'[-t int]  : temporal indice (default: 1)'
130      print*,'[-z int]  : vertical indice (default: 1)'
131      print*,'[-mt int] : extent of temporal mean (default: 0)'
132      print*,'[-mz int] : extent of vertical mean (default: 0)'
133      print*,'[-o str]  : output spectra file name (default: spectra)'
134      print*,'[-u str]  : name of zonal wind in netcdf file'
135      print*,'[-v str]  : name of meridional wind in netcdf file'
136      print*,'[-lat str]        : name of latitude field in netcdf file'
137      print*,'[-lon str]        : name of longitude field in netcdf file'
138      print*,'[-alt str]        : name of altitude field in netcdf file. ''none'' if no altitude axis.'
139      print*,'[-time str]       : name of time field in netcdf file. ''none'' if no time axis.'
140      print*,'[-divrot] : total, divergence and vorticity spectra coefficient in output file'
141      print*,''
142      print*,'Compute the kinetic energy spectrum of a wind field on a longitude-latitude grid.'
143      print*,'The grid must have a redundant point in longitude.'
144      print*,'Ideally the analysis works better on a symetric grid (2N points in longitude and N points in latitude).'
145      print*,'The output text file (called spectra by default) give the decomposition'
146      print*,'of the velocity on the vectorial spherical harmonic basis.'
147      stop 'End help'
148    else
149      file_netcdf = arg(i)
150      i = i + 1
151  end if
152end do
153
154!**********
155!Check input/output files
156inquire(file=trim(file_netcdf),exist=is_file)
157if (.not. is_file) then
158  print*,"no netcdf file: ",trim(file_netcdf)
159  stop "Stopped"
160end if
161call system('rm -f '//trim(file_spectra))
162
163!**********
164!Netcdf file dimensions
165ierror=nf90_open(trim(file_netcdf),NF90_NOWRITE,idfile)
166
167ierror=-99999
168i=0
169ierror=nf90_inq_dimid(idfile,trim(name_lat),idlat) !try user named latitude
170do while (ierror /= 0 .and. i <= size(tab_name_lat)-1)!try automatic named latitude
171  i=i+1
172  ierror=nf90_inq_dimid(idfile,trim(tab_name_lat(i)),idlat)
173end do
174if (ierror == 0) then
175  ierror=nf90_inquire_dimension(idfile,idlat,tmp_char,nlat)
176else
177  print*,'latitude dimension not found!'
178  stop ' must have this... use ''-lat latitudefieldname'' option.'
179end if
180
181ierror=-99999
182i=0
183ierror=nf90_inq_dimid(idfile,trim(name_long),idlong) !try user named longitude
184do while (ierror /= 0 .and. i <= size(tab_name_long)-1)  !try automatic named longitude
185  i=i+1
186  ierror=nf90_inq_dimid(idfile,trim(tab_name_long(i)),idlong)
187end do
188if (ierror == 0) then
189  ierror=nf90_inquire_dimension(idfile,idlong,tmp_char,nlongp1)
190else
191  print*,'longitude dimension not found!'
192  stop ' must have this... use ''-lon longitudefieldname'' option.'
193end if
194
195ierror=-99999
196i=0
197ierror=nf90_inq_dimid(idfile,trim(name_alt),idalt) !try user named altitude
198do while (ierror /= 0 .and. i <= size(tab_name_alt)-1)!try automatic named altitude
199  i=i+1
200  ierror=nf90_inq_dimid(idfile,trim(tab_name_alt(i)),idalt)
201end do
202if (ierror == 0) then
203  ierror=nf90_inquire_dimension(idfile,idalt,tmp_char,nalt)
204else
205  if (name_alt == 'none') then
206    print*,'netcdf file without altitude axis: all the same.'
207  else
208    print*,'altitude dimension not found!'
209    stop ' if no altitude axis, use ''-alt none'' option...'
210  end if
211end if
212
213ierror=-99999
214i=0
215ierror=nf90_inq_dimid(idfile,trim(name_time),idtime) !try user named time
216do while (ierror /= 0 .and. i <= size(tab_name_time)-1)  !try automatic named time
217  i=i+1
218  ierror=nf90_inq_dimid(idfile,trim(tab_name_time(i)),idtime)
219end do
220if (ierror == 0) then
221  ierror=nf90_inquire_dimension(idfile,idtime,tmp_char,ntime)
222else
223  if (name_time == 'none') then
224    print*,'netcdf file without time axis: all the same.'
225  else
226    print*,'time dimension not found!'
227    stop ' if no time axis, use ''-time none'' option...'
228  end if
229end if
230
231!**********
232!Check input indices and add average indices
233if (n_indt == 0) n_indt = 1
234tmp_tab_int(:) = tab_indt(:)
235deallocate(tab_indt)
236allocate(tab_indt(n_indt*(meant+1)))
237do t = 1,n_indt
238do mt = 0,meant
239  tab_indt((t-1)*(meant+1)+1+mt) = tmp_tab_int(t)+mt
240end do
241end do 
242do t = 1,n_indt*(meant+1)
243  if (tab_indt(t) > ntime .or. tab_indt(t) < 1) tab_indt(t) = 1
244end do
245
246if (n_indz == 0) n_indz = 1
247tmp_tab_int(:) = tab_indz(:)
248deallocate(tab_indz)
249allocate(tab_indz(n_indz*(meanz+1)))
250do z = 1,n_indz
251do mz = 0,meanz
252  tab_indz((z-1)*(meanz+1)+1+mz) = tmp_tab_int(z)+mz
253end do
254end do 
255do z = 1,n_indz*(meanz+1)
256  if (tab_indz(z) > nalt .or. tab_indz(z) < 1) tab_indz(z) = 1
257end do
258
259!**********
260!Loop over time and altitude indices
261allocate(spectra_config(nlat,n_indt*n_indz))
262allocate(spectra_b_config(nlat,n_indt*n_indz))
263allocate(spectra_c_config(nlat,n_indt*n_indz))
264allocate(spectra_global(nlat,n_indt*(meant+1)*n_indz*(meanz+1)))
265allocate(spectra_b_global(nlat,n_indt*(meant+1)*n_indz*(meanz+1)))
266allocate(spectra_c_global(nlat,n_indt*(meant+1)*n_indz*(meanz+1)))
267do t = 1,n_indt
268do z = 1,n_indz
269number_config = (t-1)*n_indz+z
270do mt = 0,meant
271do mz = 0,meanz
272number_local = mt*(meanz+1)+mz+1
273number_global = (number_config-1)*(meant+1)*(meanz+1) + number_local
274if (first_reading) then
275  print*,'First netcdf file reading...'
276  first_reading = .false.
277end if
278
279!**********
280!Netcdf file reading
281ierror=-99999
282i=0
283ierror=nf90_inq_varid(idfile,trim(name_u),idzonal_wind) !try user named zonal wind
284do while (ierror /= 0 .and. i <= size(tab_name_u)-1)        !try automatic named zonal wind
285  i=i+1
286  ierror=nf90_inq_varid(idfile,trim(tab_name_u(i)),idzonal_wind)
287end do
288if (ierror == 0) then
289  allocate(zonal_wind_lmdz(nlongp1,nlat))
290  if (name_alt == 'none') then
291    if(name_time == 'none') then
292      ierror=nf90_get_var(idfile,idzonal_wind,zonal_wind_lmdz,(/1,1/),(/nlongp1,nlat/))
293    else
294      ierror=nf90_get_var(idfile,idzonal_wind,zonal_wind_lmdz,(/1,1,tab_indt((t-1)*(meant+1)+mt+1)/),(/nlongp1,nlat,1/))
295    end if
296  end if
297  if (name_alt /= 'none') then
298    if(name_time == 'none') then
299      ierror=nf90_get_var(idfile,idzonal_wind,zonal_wind_lmdz,(/1,1,tab_indz((z-1)*(meanz+1)+1)/),(/nlongp1,nlat,1/))
300    else
301      ierror=nf90_get_var(idfile,idzonal_wind,zonal_wind_lmdz,&
302                                (/1,1,tab_indz((z-1)*(meanz+1)+1),tab_indt((t-1)*(meant+1)+mt+1)/),(/nlongp1,nlat,1,1/))
303    end if
304  end if
305else
306  print*,'zonal wind variable not found!'
307  stop ' must have this... use ''-u zonalwindfieldname'' option.'
308end if
309
310ierror=-99999
311i=0
312ierror=nf90_inq_varid(idfile,trim(name_v),idmerid_wind) !try user named meridional wind
313do while (ierror /= 0 .and. i <= size(tab_name_v)-1)        !try automatic named meridional wind
314  i=i+1
315  ierror=nf90_inq_varid(idfile,trim(tab_name_v(i)),idmerid_wind)
316end do
317if (ierror == 0) then
318  allocate(merid_wind_lmdz(nlongp1,nlat))
319  if (name_alt == 'none') then
320    if(name_time == 'none') then
321      ierror=nf90_get_var(idfile,idmerid_wind,merid_wind_lmdz,(/1,1/),(/nlongp1,nlat/))
322    else
323      ierror=nf90_get_var(idfile,idmerid_wind,merid_wind_lmdz,(/1,1,tab_indt((t-1)*(meant+1)+mt+1)/),(/nlongp1,nlat,1/))
324    end if
325  end if
326  if (name_alt /= 'none') then
327    if(name_time == 'none') then
328      ierror=nf90_get_var(idfile,idmerid_wind,merid_wind_lmdz,(/1,1,tab_indz((z-1)*(meanz+1)+1)/),(/nlongp1,nlat,1/))
329    else
330      ierror=nf90_get_var(idfile,idmerid_wind,merid_wind_lmdz,&
331                                (/1,1,tab_indz((z-1)*(meanz+1)+1),tab_indt((t-1)*(meant+1)+mt+1)/),(/nlongp1,nlat,1,1/))
332    end if
333  end if
334else
335  print*,'meridional wind variable not found!'
336  stop ' must have this... use ''-v meridionalwindfieldname'' option.'
337end if
338
339!**********From lmdz to spherepack convention
340nlong=nlongp1-1
341allocate(zonal_wind(nlat,nlong))
342allocate(merid_wind(nlat,nlong))
343zonal_wind(:,:)=transpose(zonal_wind_lmdz(1:nlong,:))
344merid_wind(:,:)=transpose(merid_wind_lmdz(1:nlong,:))
345!if (mod(nlat,2) == 0) then
346!  merid_wind(1:nlat/2,:)=-merid_wind(1:nlat/2,:)
347!else
348!  merid_wind(1:(nlat+1)/2,:)=-merid_wind(1:(nlat+1)/2,:)
349!end if
350
351!#####
352!Spectra computation
353!#####
354!**********Maximum value for m
355if (mod(nlong,2) == 0) then
356  mmax = min(nlat,nlong/2)
357else
358  mmax = min(nlat,(nlong+1)/2)
359end if
360
361!**********Vhaeci function (initialisations for Vhaec function)
362if (mod(nlong,2) == 0) then
363  l1 = min(nlat,nlong/2)
364else
365  l1 = min(nlat,(nlong+1)/2)
366end if
367if (mod(nlat,2) == 0) then
368  l2 = nlat/2
369else
370  l2 = (nlat+1)/2
371end if
372lvhaec=4*nlat*l2+3*max(l1-2,0)*(nlat+nlat-l1-1)+nlong+15
373allocate(wvhaec(lvhaec))
374wvhaec(:) = 0.
375ldwork=2*(nlat+2)
376allocate(dwork(ldwork))
377dwork(:) = 0.
378ierror=3
379call vhaeci(nlat,nlong,wvhaec,lvhaec,dwork,ldwork,ierror)
380 
381!**********Vhaeci function result
382select case (ierror)
383  case(1)
384    print*,'Vhaeci: ERROR on nlat'
385  case(2)
386    print*,'Vhaeci: ERROR on nlong'
387  case(3)
388    print*,'Vhaeci: ERROR on lvhaec'
389  case(4)
390    print*,'Vhaeci: ERROR on ldwork'
391end select
392
393!**********Vhaec function (computes spectra coefficients)
394idvw=nlat
395jdvw=nlong
396mdab=mmax
397ndab=nlat
398allocate(br(mdab,ndab))
399allocate(bi(mdab,ndab))
400allocate(cr(mdab,ndab))
401allocate(ci(mdab,ndab))
402br(:,:) = 0.
403bi(:,:) = 0.
404cr(:,:) = 0.
405ci(:,:) = 0.
406if (mod(nlong,2) == 0) then
407  l1 = min(nlat,nlong/2)
408else
409  l1 = min(nlat,(nlong+1)/2)
410end if
411if (mod(nlat,2) == 0) then
412  l2 = nlat/2
413else
414  l2 = (nlat+1)/2
415end if
416lvhaec=4*nlat*l2+3*max(l1-2,0)*(nlat+nlat-l1-1)+nlong+15
417if (ityp .le. 2) then
418  lwork=nlat*(2*nt*nlong+max(6*l2,nlong))
419else
420  lwork=l2*(2*nt*nlong+max(6*nlat,nlong))
421end if
422allocate(work(lwork))
423work(:) = 0.
424ierror=3
425call vhaec(nlat,nlong,ityp,nt,merid_wind,zonal_wind,idvw,jdvw,br,bi,cr,ci,mdab,ndab,wvhaec,lvhaec,work,lwork,ierror)
426
427!**********Vhaec function result
428select case (ierror)
429  case(1)
430    print*,'Vhaec: ERROR on nlat'
431  case(2)
432    print*,'Vhaec: ERROR on nlong'
433  case(3)
434    print*,'Vhaec: ERROR on ityp'
435  case(4)
436    print*,'Vhaec: ERROR on nt'
437  case(5)
438    print*,'Vhaec: ERROR on idvw'
439  case(6)
440    print*,'Vhaec: ERROR on jdvw'
441  case(7)
442    print*,'Vhaec: ERROR on mdab'
443  case(8)
444    print*,'Vhaec: ERROR on ndab'
445  case(9)
446    print*,'Vhaec: ERROR on lvhaec'
447  case(10)
448    print*,'Vhaec: ERROR on lwork'
449end select
450
451!**********Energy spectra
452allocate(norm(mMax,nLat))
453allocate(norm_b(mMax,nLat))
454allocate(norm_c(mMax,nLat))
455do i = 1,mMax
456do j = 1,nLat
457  norm(i,j)=(br(i,j)*br(i,j)+bi(i,j)*bi(i,j)+cr(i,j)*cr(i,j)+ci(i,j)*ci(i,j))
458  norm_b(i,j)=(br(i,j)*br(i,j)+bi(i,j)*bi(i,j))
459  norm_c(i,j)=(cr(i,j)*cr(i,j)+ci(i,j)*ci(i,j))
460end do
461end do
462allocate(spectra(nLat))
463allocate(spectra_b(nLat))
464allocate(spectra_c(nLat))
465do j = 1,nLat
466  spectra(j) = 0.5*norm(1,j)
467  spectra_b(j) = 0.5*norm_b(1,j)
468  spectra_c(j) = 0.5*norm_c(1,j)
469  do i = 2,mMax
470    spectra(j) = spectra(j) + norm(i,j)
471    spectra_b(j) = spectra_b(j) + norm_b(i,j)
472    spectra_c(j) = spectra_c(j) + norm_c(i,j)
473  end do
474end do
475!do j = 2,nLat
476!  spectra(j) = spectra(j)/((j-1)*j)
477!  spectra_b(j) = spectra_b(j)/((j-1)*j)
478!  spectra_c(j) = spectra_c(j)/((j-1)*j)
479!end do
480write(*,'(a16,i5,a3,i5,a1)') 'done spectra (t=',tab_indt((t-1)*(meant+1)+mt+1),',z=',tab_indz((z-1)*(meanz+1)+mz+1),')'
481
482!**********Storage of all spectra
483spectra_global(:,number_global) = spectra(:)
484spectra_b_global(:,number_global) = spectra_b(:)
485spectra_c_global(:,number_global) = spectra_c(:)
486
487!**********Some cleaning
488deallocate(br)
489deallocate(bi)
490deallocate(cr)
491deallocate(ci)
492deallocate(wvhaec)
493deallocate(dwork)
494deallocate(work)
495deallocate(norm)
496deallocate(norm_b)
497deallocate(norm_c)
498deallocate(spectra)
499deallocate(spectra_b)
500deallocate(spectra_c)
501deallocate(zonal_wind)
502deallocate(merid_wind)
503deallocate(zonal_wind_lmdz)
504deallocate(merid_wind_lmdz)
505
506!**********
507!End of loop over time and altitude indices
508end do
509end do
510write(*,'(a24,i5,a3,i5,a8)') 'spectra computed for (t=',tab_indt((t-1)*(meant+1)+1),',z=',tab_indz((z-1)*(meanz+1)+1),') config'
511write(*,*) '**********'
512end do
513end do
514
515!**********
516!Spectra mean
517do t = 1,n_indt
518do z = 1,n_indz
519number_config = (t-1)*n_indz+z
520spectra_config(:,number_config) = 0.
521do mt = 0,meant
522do mz = 0,meanz
523number_local = mt*(meanz+1)+mz+1
524number_global = (number_config-1)*(meant+1)*(meanz+1) + number_local
525spectra_config(:,number_config) = spectra_config(:,number_config) + &
526                spectra_global(:,number_global)
527spectra_b_config(:,number_config) = spectra_b_config(:,number_config) + &
528                spectra_b_global(:,number_global)
529spectra_c_config(:,number_config) = spectra_c_config(:,number_config) + &
530                spectra_c_global(:,number_global)
531end do
532end do
533spectra_config(:,number_config) = spectra_config(:,number_config)/(meanz+1)/(meant+1)
534spectra_b_config(:,number_config) = spectra_b_config(:,number_config)/(meanz+1)/(meant+1)
535spectra_c_config(:,number_config) = spectra_c_config(:,number_config)/(meanz+1)/(meant+1)
536end do
537end do
538
539!**********Writing header of file_spectra
540open(unit=10,file=file_spectra,action="write",position="rewind")
541write(10,'(a10,a2)',advance='no') '#Spherical',' '
542do t = 1,n_indt
543do z = 1,n_indz
544  write(10,'(a50)',advance='no') file_netcdf
545end do
546end do
547write(10,*)
548write(10,'(a11,a1)',advance='no') '#wavenumber',' '
549do t = 1,n_indt
550do z = 1,n_indz
551  write(10,'(a2,i5,a3,i5,a35)',advance='no') 't=',tab_indt((t-1)*(meant+1)+1),' z=',tab_indz((z-1)*(meanz+1)+1),' '
552end do
553end do
554write(10,*)
555write(10,'(a1,a11)',advance='no') '#',' '
556do t = 1,n_indt
557do z = 1,n_indz
558if (is_div_rot) then
559  write(10,'(a8,a8,a8,a8,a8,a10)',advance='no') 'spec_tot',' ','spec_div',' ','spec_rot',' '
560else
561  write(10,'(a8,a42)',advance='no') 'spec_tot',' '
562end if
563end do
564end do
565write(10,*)
566!**********Writing file_spectra
567do j=1,nLat
568  write(10,'(i5,a6)',advance='no') j-1,' '
569  do t = 1,n_indt
570  do z = 1,n_indz
571  number_config = (t-1)*n_indz+z
572  if (is_div_rot) then
573    write(10,'(e13.6E2,a3,e13.6E2,a3,e13.6E2,a5)',advance='no') spectra_config(j,number_config),' ',&
574                        spectra_b_config(j,number_config),' ',spectra_c_config(j,number_config),' '
575  else
576    write(10,'(e13.6E2,a37)',advance='no') spectra_config(j,number_config),' '
577  end if
578  end do
579  end do
580  if (j /= nlat) write(10,*)
581end do
582close(10)
583
584print *, "***SUCCESS writing ",trim(file_spectra)
585
586!***********Some cleaning
587deallocate(spectra_config)
588deallocate(spectra_b_config)
589deallocate(spectra_c_config)
590deallocate(spectra_global)
591deallocate(spectra_b_global)
592deallocate(spectra_c_global)
593deallocate(tab_indt)
594deallocate(tab_indz)
595deallocate(tmp_tab_int)
596deallocate(arg)
597
598contains
599  subroutine check(status)
600    integer, intent (in) :: status
601   
602    if(status /= nf90_noerr) then
603      print *, trim(nf90_strerror(status))
604      stop "Stopped"
605    end if
606  end subroutine check 
607
608END PROGRAM spectra_analysis
Note: See TracBrowser for help on using the repository browser.