[816] | 1 | program fft |
---|
| 2 | |
---|
| 3 | ! SL 01/2010: |
---|
| 4 | ! This program reads 4D (lon-lat-alt-time) fields recast in log P coordinates |
---|
| 5 | ! |
---|
| 6 | ! it computes fft of temperature, zonal and merid winds from high-frequency outputs: |
---|
| 7 | ! |
---|
| 8 | ! fftaT -- 4D -- FFT in amplitude of temperature field (K) |
---|
| 9 | ! fftau -- 4D -- FFT in amplitude of zonal wind (m s-1) |
---|
| 10 | ! fftav -- 4D -- FFT in amplitude of meridional wind (m s-1) |
---|
| 11 | ! ulf -- 4D -- low freq part of zonal wind perturbation uprim (m s-1) |
---|
| 12 | ! ubf -- 4D -- band freq part of zonal wind perturbation uprim (m s-1) |
---|
| 13 | ! uhf -- 4D -- high freq part of zonal wind perturbation uprim (m s-1) |
---|
| 14 | ! vlf -- 4D -- low freq part of meridional wind perturbation vprim (m s-1) |
---|
| 15 | ! vbf -- 4D -- band freq part of meridional wind perturbation vprim (m s-1) |
---|
| 16 | ! vhf -- 4D -- high freq part of meridional wind perturbation vprim (m s-1) |
---|
| 17 | ! wlf -- 4D -- low freq part of vertical wind perturbation wprim (Pa s-1) |
---|
| 18 | ! wbf -- 4D -- band freq part of vertical wind perturbation wprim (Pa s-1) |
---|
| 19 | ! whf -- 4D -- high freq part of vertical wind perturbation wprim (Pa s-1) |
---|
| 20 | ! Tlf -- 4D -- low freq part of temperature perturbation Tprim (K) |
---|
| 21 | ! Tbf -- 4D -- band freq part of temperature perturbation Tprim (K) |
---|
| 22 | ! Thf -- 4D -- high freq part of temperature perturbation Tprim (K) |
---|
| 23 | ! |
---|
| 24 | ! Minimal requirements and dependencies: |
---|
| 25 | ! The dataset must include the following data: |
---|
| 26 | ! - pressure vertical coordinate |
---|
| 27 | ! - atmospheric temperature |
---|
| 28 | ! - zonal, meridional and vertical winds |
---|
| 29 | ! |
---|
| 30 | ! We use the FFTW library: http://www.fftw.org |
---|
| 31 | ! These routines are in C, but also include Fortran interfaces. |
---|
| 32 | ! |
---|
| 33 | ! Convention: qbar <=> zonal average / qstar = q - qbar |
---|
| 34 | ! qmean <=> temporal average / qprim = q - qmean |
---|
| 35 | ! |
---|
| 36 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 37 | ! FILTRES |
---|
| 38 | !!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 39 | ! low frequencies: qlf= lowfreq(qprim) |
---|
| 40 | ! band frequencies: qbf=bandfreq(qprim) |
---|
| 41 | ! high frequencies: qhf=highfreq(qprim) |
---|
| 42 | ! |
---|
| 43 | ! Les frequences seuils sont ajustables dans filter.h |
---|
| 44 | ! |
---|
| 45 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 46 | |
---|
| 47 | implicit none |
---|
| 48 | |
---|
| 49 | include "netcdf.inc" ! NetCDF definitions |
---|
| 50 | |
---|
| 51 | character (len=128) :: infile ! input file name (name_P.nc) |
---|
| 52 | character (len=128) :: outfile1,outfile2,outfile3,outfile4 ! output file names |
---|
| 53 | |
---|
| 54 | character (len=64) :: text ! to store some text |
---|
| 55 | integer infid ! NetCDF input file ID |
---|
| 56 | integer outfid1,outfid2,outfid3,outfid4 ! NetCDF output files ID |
---|
| 57 | integer lon_dimid1,lat_dimid1,alt_dimid1,time_dimid1 ! NetCDF dimension IDs |
---|
| 58 | integer lon_dimid2,lat_dimid2,alt_dimid2,time_dimid2 ! NetCDF dimension IDs |
---|
| 59 | integer lon_dimid3,lat_dimid3,alt_dimid3,time_dimid3 ! NetCDF dimension IDs |
---|
| 60 | integer lon_dimid4,lat_dimid4,alt_dimid4,time_dimid4 ! NetCDF dimension IDs |
---|
| 61 | integer lon_varid,lat_varid,alt_varid,time_varid |
---|
| 62 | integer :: datashape1d ! shape of 1D datasets |
---|
| 63 | integer,dimension(4) :: datashape4d ! shape of 4D datasets |
---|
| 64 | |
---|
| 65 | real :: miss_val=9.99e+33 ! special "missing value" to specify missing data |
---|
| 66 | real,parameter :: miss_val_def=9.99e+33 ! default value for "missing value" |
---|
| 67 | real :: pi |
---|
| 68 | real,dimension(:),allocatable :: lon ! longitude |
---|
| 69 | integer lonlength ! # of grid points along longitude |
---|
| 70 | real,dimension(:),allocatable :: lat ! latitude |
---|
| 71 | real,dimension(:),allocatable :: latrad ! latitude in radian |
---|
| 72 | integer latlength ! # of grid points along latitude |
---|
| 73 | real,dimension(:),allocatable :: plev ! Pressure levels (Pa) |
---|
| 74 | integer altlength ! # of grid point along altitude (of input datasets) |
---|
| 75 | real,dimension(:),allocatable :: time ! time |
---|
| 76 | real,dimension(:),allocatable :: freq ! frequencies of the FFT (only timelength/2+1 values) |
---|
| 77 | integer timelength ! # of points along time |
---|
| 78 | real,dimension(:,:,:,:),allocatable :: temp ! atmospheric temperature |
---|
| 79 | real,dimension(:,:,:,:),allocatable :: vitu ! zonal wind (in m/s) |
---|
| 80 | real,dimension(:,:,:,:),allocatable :: vitv ! meridional wind (in m/s) |
---|
| 81 | real,dimension(:,:,:,:),allocatable :: vitw ! vertical wind (in Pa/s) |
---|
| 82 | |
---|
| 83 | !!! output variables |
---|
| 84 | real,dimension(:,:,:,:),allocatable :: fftaT ! FFT in amplitude of temperature (K) |
---|
| 85 | real,dimension(:,:,:,:),allocatable :: fftau ! FFT in amplitude of zonal wind (m s-1) |
---|
| 86 | real,dimension(:,:,:,:),allocatable :: fftav ! FFT in amplitude of meridional wind (m s-1) |
---|
| 87 | real,dimension(:,:,:,:),allocatable :: fftaw ! FFT in amplitude of vertical wind (Pa s-1) |
---|
| 88 | real,dimension(:,:,:,:),allocatable :: ulf ! low freq part of zonal wind perturbation uprim (m s-1) |
---|
| 89 | real,dimension(:,:,:,:),allocatable :: ubf ! band freq part of zonal wind perturbation uprim (m s-1) |
---|
| 90 | real,dimension(:,:,:,:),allocatable :: uhf ! high freq part of zonal wind perturbation uprim (m s-1) |
---|
| 91 | real,dimension(:,:,:,:),allocatable :: vlf ! low freq part of meridional wind perturbation vprim (m s-1) |
---|
| 92 | real,dimension(:,:,:,:),allocatable :: vbf ! band freq part of meridional wind perturbation vprim (m s-1) |
---|
| 93 | real,dimension(:,:,:,:),allocatable :: vhf ! high freq part of meridional wind perturbation vprim (m s-1) |
---|
| 94 | real,dimension(:,:,:,:),allocatable :: wlf ! low freq part of vertical wind perturbation vprim (Pa s-1) |
---|
| 95 | real,dimension(:,:,:,:),allocatable :: wbf ! band freq part of vertical wind perturbation vprim (Pa s-1) |
---|
| 96 | real,dimension(:,:,:,:),allocatable :: whf ! high freq part of vertical wind perturbation vprim (Pa s-1) |
---|
| 97 | real,dimension(:,:,:,:),allocatable :: Tlf ! low freq part of temperature perturbation Tprim (K) |
---|
| 98 | real,dimension(:,:,:,:),allocatable :: Tbf ! band freq part of temperature perturbation Tprim (K) |
---|
| 99 | real,dimension(:,:,:,:),allocatable :: Thf ! high freq part of temperature perturbation Tprim (K) |
---|
| 100 | |
---|
| 101 | ! local variables |
---|
| 102 | real,dimension(:,:,:,:),allocatable :: uprim |
---|
| 103 | real,dimension(:,:,:,:),allocatable :: vprim |
---|
| 104 | real,dimension(:,:,:,:),allocatable :: wprim |
---|
| 105 | real,dimension(:,:,:,:),allocatable :: Tprim |
---|
| 106 | |
---|
| 107 | ! lon,lat,alt |
---|
| 108 | real,dimension(:,:,:),allocatable :: umean |
---|
| 109 | real,dimension(:,:,:),allocatable :: vmean |
---|
| 110 | real,dimension(:,:,:),allocatable :: wmean |
---|
| 111 | real,dimension(:,:,:),allocatable :: Tmean |
---|
| 112 | |
---|
| 113 | ! for FFTW routines |
---|
| 114 | real,dimension(:),allocatable :: wndow |
---|
| 115 | double precision,dimension(:),allocatable :: var,fltvar |
---|
| 116 | double complex,dimension(:),allocatable :: fftvar,fltfft |
---|
| 117 | double complex,dimension(:),allocatable :: filtrelf,filtrebf,filtrehf |
---|
| 118 | integer :: M_fft |
---|
| 119 | integer*8 :: planf,planb |
---|
| 120 | |
---|
| 121 | |
---|
| 122 | integer ierr,ierr1,ierr2 ! NetCDF routines return codes |
---|
| 123 | integer i,j,ilon,ilat,ilev,itim ! for loops |
---|
| 124 | logical flagfft |
---|
| 125 | logical :: lmdflag |
---|
| 126 | |
---|
| 127 | include "planet.h" |
---|
| 128 | include "filter.h" |
---|
| 129 | |
---|
| 130 | #include <fftw3.f> |
---|
| 131 | |
---|
| 132 | !=============================================================================== |
---|
| 133 | ! 1. Input parameters |
---|
| 134 | !=============================================================================== |
---|
| 135 | |
---|
| 136 | pi = 2.*asin(1.) |
---|
| 137 | |
---|
| 138 | write(*,*) "" |
---|
| 139 | write(*,*) "You are working on the atmosphere of ",planet |
---|
| 140 | |
---|
| 141 | !=============================================================================== |
---|
| 142 | ! 1.1 Input file |
---|
| 143 | !=============================================================================== |
---|
| 144 | |
---|
| 145 | write(*,*) "" |
---|
| 146 | write(*,*) " Program valid for files with pressure axis (*_P.nc)" |
---|
| 147 | write(*,*) "Enter input file name:" |
---|
| 148 | |
---|
| 149 | read(*,'(a128)') infile |
---|
| 150 | write(*,*) "" |
---|
| 151 | |
---|
| 152 | ! open input file |
---|
| 153 | |
---|
| 154 | ierr = NF_OPEN(infile,NF_NOWRITE,infid) |
---|
| 155 | if (ierr.ne.NF_NOERR) then |
---|
| 156 | write(*,*) 'ERROR: Pb opening file ',trim(infile) |
---|
| 157 | stop "" |
---|
| 158 | endif |
---|
| 159 | |
---|
| 160 | !=============================================================================== |
---|
| 161 | ! 1.2 Get grids in lon,lat,alt(pressure),time |
---|
| 162 | !=============================================================================== |
---|
| 163 | |
---|
| 164 | call get_iddim(infid,lat_varid,latlength,lon_varid,lonlength,& |
---|
| 165 | alt_varid,altlength,time_varid,timelength,lmdflag ) |
---|
| 166 | |
---|
| 167 | allocate(lon(lonlength)) |
---|
| 168 | ierr=NF_GET_VAR_REAL(infid,lon_varid,lon) |
---|
| 169 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading longitude" |
---|
| 170 | |
---|
| 171 | allocate(lat(latlength)) |
---|
| 172 | ierr=NF_GET_VAR_REAL(infid,lat_varid,lat) |
---|
| 173 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading lat" |
---|
| 174 | |
---|
| 175 | allocate(latrad(latlength)) |
---|
| 176 | latrad = lat*pi/180. |
---|
| 177 | |
---|
| 178 | allocate(plev(altlength)) |
---|
| 179 | ierr=NF_GET_VAR_REAL(infid,alt_varid,plev) |
---|
| 180 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading altitude (ie pressure levels)" |
---|
| 181 | |
---|
| 182 | allocate(time(timelength)) |
---|
| 183 | ierr=NF_GET_VAR_REAL(infid,time_varid,time) |
---|
| 184 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading time" |
---|
| 185 | |
---|
| 186 | !=============================================================================== |
---|
| 187 | ! 1.3 Get output file name |
---|
| 188 | !=============================================================================== |
---|
| 189 | write(*,*) "" |
---|
| 190 | !write(*,*) "Enter output file name" |
---|
| 191 | !read(*,*) outfile |
---|
| 192 | outfile1=infile(1:len_trim(infile)-3)//"_UFFT.nc" |
---|
| 193 | outfile2=infile(1:len_trim(infile)-3)//"_VFFT.nc" |
---|
| 194 | outfile3=infile(1:len_trim(infile)-3)//"_WFFT.nc" |
---|
| 195 | outfile4=infile(1:len_trim(infile)-3)//"_TFFT.nc" |
---|
| 196 | write(*,*) "Output file names are: " |
---|
| 197 | if (ok_out(1)) write(*,*) trim(outfile1) |
---|
| 198 | if (ok_out(2)) write(*,*) trim(outfile2) |
---|
| 199 | if (ok_out(3)) write(*,*) trim(outfile3) |
---|
| 200 | if (ok_out(4)) write(*,*) trim(outfile4) |
---|
| 201 | |
---|
| 202 | |
---|
| 203 | !=============================================================================== |
---|
| 204 | ! 2.1 Store needed fields |
---|
| 205 | !=============================================================================== |
---|
| 206 | |
---|
| 207 | !=============================================================================== |
---|
| 208 | ! 2.1.1 Atmospheric temperature |
---|
| 209 | !=============================================================================== |
---|
| 210 | allocate(temp(lonlength,latlength,altlength,timelength)) |
---|
| 211 | |
---|
| 212 | text="temp" |
---|
| 213 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,temp,ierr1,ierr2) |
---|
| 214 | if (ierr1.ne.NF_NOERR) then |
---|
| 215 | write(*,*) " looking for t instead... " |
---|
| 216 | text="t" |
---|
| 217 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,temp,ierr1,ierr2) |
---|
| 218 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get temperature ID" |
---|
| 219 | endif |
---|
| 220 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading temperature" |
---|
| 221 | |
---|
| 222 | !=============================================================================== |
---|
| 223 | ! 2.1.2 Winds |
---|
| 224 | !=============================================================================== |
---|
| 225 | allocate(vitu(lonlength,latlength,altlength,timelength)) |
---|
| 226 | allocate(vitv(lonlength,latlength,altlength,timelength)) |
---|
| 227 | allocate(vitw(lonlength,latlength,altlength,timelength)) |
---|
| 228 | |
---|
| 229 | ! zonal wind vitu (in m/s) |
---|
| 230 | text="vitu" |
---|
| 231 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitu,ierr1,ierr2) |
---|
| 232 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitu ID" |
---|
| 233 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zonal wind" |
---|
| 234 | |
---|
| 235 | ! meridional wind vitv (in m/s) |
---|
| 236 | text="vitv" |
---|
| 237 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitv,ierr1,ierr2) |
---|
| 238 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitv ID" |
---|
| 239 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading meridional wind" |
---|
| 240 | |
---|
| 241 | ! vertical wind vitw (in Pa/s) |
---|
| 242 | text="vitw" |
---|
| 243 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitw,ierr1,ierr2) |
---|
| 244 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitw ID" |
---|
| 245 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading vertical wind" |
---|
| 246 | |
---|
| 247 | !=============================================================================== |
---|
| 248 | ! 2.2 Computations |
---|
| 249 | !=============================================================================== |
---|
| 250 | |
---|
| 251 | print*,"debut calcul" |
---|
| 252 | |
---|
| 253 | !=============================================================================== |
---|
| 254 | ! 2.2.1 FFT and filtering |
---|
| 255 | !=============================================================================== |
---|
| 256 | ! allocations |
---|
| 257 | !------------- |
---|
| 258 | allocate(fftau(lonlength,latlength,altlength,timelength)) |
---|
| 259 | allocate(fftav(lonlength,latlength,altlength,timelength)) |
---|
| 260 | allocate(fftaw(lonlength,latlength,altlength,timelength)) |
---|
| 261 | allocate(fftaT(lonlength,latlength,altlength,timelength)) |
---|
| 262 | allocate(uprim(lonlength,latlength,altlength,timelength)) |
---|
| 263 | allocate(vprim(lonlength,latlength,altlength,timelength)) |
---|
| 264 | allocate(wprim(lonlength,latlength,altlength,timelength)) |
---|
| 265 | allocate(Tprim(lonlength,latlength,altlength,timelength)) |
---|
| 266 | allocate(ulf(lonlength,latlength,altlength,timelength)) |
---|
| 267 | allocate(vlf(lonlength,latlength,altlength,timelength)) |
---|
| 268 | allocate(wlf(lonlength,latlength,altlength,timelength)) |
---|
| 269 | allocate(Tlf(lonlength,latlength,altlength,timelength)) |
---|
| 270 | allocate(ubf(lonlength,latlength,altlength,timelength)) |
---|
| 271 | allocate(vbf(lonlength,latlength,altlength,timelength)) |
---|
| 272 | allocate(wbf(lonlength,latlength,altlength,timelength)) |
---|
| 273 | allocate(Tbf(lonlength,latlength,altlength,timelength)) |
---|
| 274 | allocate(uhf(lonlength,latlength,altlength,timelength)) |
---|
| 275 | allocate(vhf(lonlength,latlength,altlength,timelength)) |
---|
| 276 | allocate(whf(lonlength,latlength,altlength,timelength)) |
---|
| 277 | allocate(Thf(lonlength,latlength,altlength,timelength)) |
---|
| 278 | |
---|
| 279 | ! lon,lat,alt |
---|
| 280 | allocate(umean(lonlength,latlength,altlength)) |
---|
| 281 | allocate(vmean(lonlength,latlength,altlength)) |
---|
| 282 | allocate(wmean(lonlength,latlength,altlength)) |
---|
| 283 | allocate(Tmean(lonlength,latlength,altlength)) |
---|
| 284 | |
---|
| 285 | ! time / frequencies |
---|
| 286 | allocate(freq(timelength)) |
---|
| 287 | allocate(wndow(timelength)) |
---|
| 288 | allocate(var(timelength)) |
---|
| 289 | allocate(fltvar(timelength)) |
---|
| 290 | M_fft = timelength/2 |
---|
| 291 | allocate(fftvar(M_fft+1)) |
---|
| 292 | allocate(fltfft(M_fft+1)) |
---|
| 293 | allocate(filtrelf(M_fft+1)) |
---|
| 294 | allocate(filtrebf(M_fft+1)) |
---|
| 295 | allocate(filtrehf(M_fft+1)) |
---|
| 296 | |
---|
| 297 | ! intermediates |
---|
| 298 | !----------------- |
---|
| 299 | |
---|
| 300 | call moytim(lonlength,latlength,altlength,timelength,miss_val,vitu,umean) |
---|
| 301 | call moytim(lonlength,latlength,altlength,timelength,miss_val,vitv,vmean) |
---|
| 302 | call moytim(lonlength,latlength,altlength,timelength,miss_val,vitw,wmean) |
---|
| 303 | call moytim(lonlength,latlength,altlength,timelength,miss_val,temp,Tmean) |
---|
| 304 | |
---|
| 305 | do ilon=1,lonlength |
---|
| 306 | do ilat=1,latlength |
---|
| 307 | do ilev=1,altlength |
---|
| 308 | do itim=1,timelength |
---|
| 309 | if ((vitu(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 310 | (umean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 311 | uprim(ilon,ilat,ilev,itim) = vitu(ilon,ilat,ilev,itim)-umean(ilon,ilat,ilev) |
---|
| 312 | else |
---|
| 313 | uprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 314 | endif |
---|
| 315 | if ((vitv(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 316 | (vmean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 317 | vprim(ilon,ilat,ilev,itim) = vitv(ilon,ilat,ilev,itim)-vmean(ilon,ilat,ilev) |
---|
| 318 | else |
---|
| 319 | vprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 320 | endif |
---|
| 321 | if ((vitw(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 322 | (wmean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 323 | wprim(ilon,ilat,ilev,itim) = vitw(ilon,ilat,ilev,itim)-wmean(ilon,ilat,ilev) |
---|
| 324 | else |
---|
| 325 | wprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 326 | endif |
---|
| 327 | if ((temp(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 328 | (Tmean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 329 | Tprim(ilon,ilat,ilev,itim) = temp(ilon,ilat,ilev,itim)-Tmean(ilon,ilat,ilev) |
---|
| 330 | else |
---|
| 331 | Tprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 332 | endif |
---|
| 333 | enddo |
---|
| 334 | enddo |
---|
| 335 | enddo |
---|
| 336 | enddo ! lonlength |
---|
| 337 | |
---|
| 338 | ! fft intermediates |
---|
| 339 | !------------- |
---|
| 340 | |
---|
| 341 | ! Define the frequencies |
---|
| 342 | do itim=1,M_fft+1 |
---|
| 343 | freq(itim) = (itim-1)/(timelength*(time(2)-time(1))) |
---|
| 344 | enddo |
---|
| 345 | do itim=M_fft+2,timelength |
---|
| 346 | freq(itim) = 0. |
---|
| 347 | enddo |
---|
| 348 | |
---|
| 349 | ! Define the window (triangle) |
---|
| 350 | do itim=1,timelength |
---|
| 351 | ! N window: |
---|
| 352 | ! wndow(itim)= 1. |
---|
| 353 | ! triangulaire de moyenne = 1 |
---|
| 354 | wndow(itim)= 2.*(1. - abs(real(itim-0.5-M_fft)/real(M_fft))) |
---|
| 355 | enddo |
---|
| 356 | |
---|
| 357 | ! Define the filters |
---|
| 358 | |
---|
| 359 | ! IN filter.h : |
---|
| 360 | ! Low cutting frequency, in Hz : fcoup1 |
---|
| 361 | ! High cutting frequency, in Hz : fcoup2 |
---|
| 362 | ! Half-width of the filters, in Hz : width |
---|
| 363 | |
---|
| 364 | !print*,"Low cutting frequency, in Hz ?" |
---|
| 365 | !read(*,*) fcoup1 |
---|
| 366 | !print*,"High cutting frequency, in Hz ?" |
---|
| 367 | !read(*,*) fcoup2 |
---|
| 368 | !print*,"Half-width of the filters, in Hz ?" |
---|
| 369 | !read(*,*) width |
---|
| 370 | !width = 3./time(timelength) |
---|
| 371 | |
---|
| 372 | do itim=1,M_fft+1 |
---|
| 373 | if (freq(itim).lt.(fcoup1-width)) then |
---|
| 374 | filtrelf(itim) = 1. |
---|
| 375 | elseif (freq(itim).gt.(fcoup1+width)) then |
---|
| 376 | filtrelf(itim) = 0. |
---|
| 377 | else |
---|
| 378 | filtrelf(itim) = (1.+sin(pi*(fcoup1-freq(itim))/(2.*width)))/2. |
---|
| 379 | endif |
---|
| 380 | if (freq(itim).lt.(fcoup2-width)) then |
---|
| 381 | filtrehf(itim) = 0. |
---|
| 382 | elseif (freq(itim).gt.(fcoup2+width)) then |
---|
| 383 | filtrehf(itim) = 1. |
---|
| 384 | else |
---|
| 385 | filtrehf(itim) = (1.-sin(pi*(fcoup2-freq(itim))/(2.*width)))/2. |
---|
| 386 | endif |
---|
| 387 | filtrebf(itim) = (1.-filtrelf(itim))*(1.-filtrehf(itim)) |
---|
| 388 | enddo |
---|
| 389 | |
---|
| 390 | |
---|
| 391 | ! fft and filtering |
---|
| 392 | !------------- |
---|
| 393 | |
---|
| 394 | !---FFTW routines |
---|
| 395 | call dfftw_plan_dft_r2c_1d(planf,timelength,var,fftvar,FFTW_MEASURE) |
---|
| 396 | call dfftw_plan_dft_c2r_1d(planb,timelength,fltfft,fltvar,FFTW_MEASURE) |
---|
| 397 | !--- |
---|
| 398 | |
---|
| 399 | do ilon=1,lonlength |
---|
| 400 | do ilat=1,latlength |
---|
| 401 | do ilev=1,altlength |
---|
| 402 | |
---|
| 403 | ! For zonal wind field |
---|
| 404 | if (ok_out(1)) then |
---|
| 405 | |
---|
| 406 | flagfft=.true. |
---|
| 407 | do itim=1,timelength |
---|
| 408 | if (uprim(ilon,ilat,ilev,itim).eq.miss_val) flagfft=.false. |
---|
| 409 | enddo |
---|
| 410 | |
---|
| 411 | if (flagfft) then |
---|
| 412 | |
---|
| 413 | ! 1/ windowing to improve spectral analysis |
---|
| 414 | var(:)=uprim(ilon,ilat,ilev,:)*wndow(:) |
---|
| 415 | ! 2/ FFT computation |
---|
| 416 | !---FFTW routines |
---|
| 417 | call dfftw_execute_dft_r2c(planf,var,fftvar) |
---|
| 418 | !--- |
---|
| 419 | ! 3/ Amplitude of the FFT, for spectral analysis |
---|
| 420 | fftau(ilon,ilat,ilev,1)=abs(fftvar(1))/M_fft |
---|
| 421 | do itim=2,M_fft |
---|
| 422 | fftau(ilon,ilat,ilev,itim) = abs(fftvar(itim))/M_fft |
---|
| 423 | enddo |
---|
| 424 | fftau(ilon,ilat,ilev,M_fft+1)=abs(fftvar(M_fft+1))/M_fft |
---|
| 425 | do itim=M_fft+2,timelength |
---|
| 426 | fftau(ilon,ilat,ilev,itim) = 0. |
---|
| 427 | enddo |
---|
| 428 | |
---|
| 429 | ! 4/ filtering the FFT in three regions |
---|
| 430 | ! filtering + normalisation (low freq) |
---|
| 431 | fltfft(:) = fftvar(:)*filtrelf(:)/timelength |
---|
| 432 | ! 5/ backward FFT for each region |
---|
| 433 | !---FFTW routines |
---|
| 434 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 435 | !--- |
---|
| 436 | ! 6/ reverse the windowing |
---|
| 437 | ulf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 438 | |
---|
| 439 | ! filtering + normalisation (band freq) |
---|
| 440 | fltfft(:) = fftvar(:)*filtrebf(:)/timelength |
---|
| 441 | !---FFTW routines |
---|
| 442 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 443 | !--- |
---|
| 444 | ubf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 445 | |
---|
| 446 | ! filtering + normalisation (high freq) |
---|
| 447 | fltfft(:) = fftvar(:)*filtrehf(:)/timelength |
---|
| 448 | !---FFTW routines |
---|
| 449 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 450 | !--- |
---|
| 451 | uhf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 452 | |
---|
| 453 | else |
---|
| 454 | fftau(ilon,ilat,ilev,itim) = miss_val |
---|
| 455 | ulf(ilon,ilat,ilev,itim) = miss_val |
---|
| 456 | ubf(ilon,ilat,ilev,itim) = miss_val |
---|
| 457 | uhf(ilon,ilat,ilev,itim) = miss_val |
---|
| 458 | endif ! flagfft |
---|
| 459 | |
---|
| 460 | endif !ok_out(1) |
---|
| 461 | |
---|
| 462 | ! For meridional wind wind field |
---|
| 463 | if (ok_out(2)) then |
---|
| 464 | |
---|
| 465 | flagfft=.true. |
---|
| 466 | do itim=1,timelength |
---|
| 467 | if (vprim(ilon,ilat,ilev,itim).eq.miss_val) flagfft=.false. |
---|
| 468 | enddo |
---|
| 469 | |
---|
| 470 | if (flagfft) then |
---|
| 471 | |
---|
| 472 | ! 1/ windowing to improve spectral analysis |
---|
| 473 | var(:)=vprim(ilon,ilat,ilev,:)*wndow(:) |
---|
| 474 | ! 2/ FFT computation |
---|
| 475 | !---FFTW routines |
---|
| 476 | call dfftw_execute_dft_r2c(planf,var,fftvar) |
---|
| 477 | !--- |
---|
| 478 | ! 3/ Amplitude of the FFT, for spectral analysis |
---|
| 479 | fftav(ilon,ilat,ilev,1)=abs(fftvar(1))/M_fft |
---|
| 480 | do itim=2,M_fft |
---|
| 481 | fftav(ilon,ilat,ilev,itim) = abs(fftvar(itim))/M_fft |
---|
| 482 | enddo |
---|
| 483 | fftav(ilon,ilat,ilev,M_fft+1)=abs(fftvar(M_fft+1))/M_fft |
---|
| 484 | do itim=M_fft+2,timelength |
---|
| 485 | fftav(ilon,ilat,ilev,itim) = 0. |
---|
| 486 | enddo |
---|
| 487 | |
---|
| 488 | ! 4/ filtering the FFT in three regions |
---|
| 489 | ! filtering + normalisation (low freq) |
---|
| 490 | fltfft(:) = fftvar(:)*filtrelf(:)/timelength |
---|
| 491 | ! 5/ backward FFT for each region |
---|
| 492 | !---FFTW routines |
---|
| 493 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 494 | !--- |
---|
| 495 | ! 6/ reverse the windowing |
---|
| 496 | vlf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 497 | |
---|
| 498 | ! filtering + normalisation (band freq) |
---|
| 499 | fltfft(:) = fftvar(:)*filtrebf(:)/timelength |
---|
| 500 | !---FFTW routines |
---|
| 501 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 502 | !--- |
---|
| 503 | vbf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 504 | |
---|
| 505 | ! filtering + normalisation (high freq) |
---|
| 506 | fltfft(:) = fftvar(:)*filtrehf(:)/timelength |
---|
| 507 | !---FFTW routines |
---|
| 508 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 509 | !--- |
---|
| 510 | vhf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 511 | |
---|
| 512 | else |
---|
| 513 | fftav(ilon,ilat,ilev,itim) = miss_val |
---|
| 514 | vlf(ilon,ilat,ilev,itim) = miss_val |
---|
| 515 | vbf(ilon,ilat,ilev,itim) = miss_val |
---|
| 516 | vhf(ilon,ilat,ilev,itim) = miss_val |
---|
| 517 | endif ! flagfft |
---|
| 518 | |
---|
| 519 | endif !ok_out(2) |
---|
| 520 | |
---|
| 521 | ! For vertical wind wind field |
---|
| 522 | if (ok_out(3)) then |
---|
| 523 | |
---|
| 524 | flagfft=.true. |
---|
| 525 | do itim=1,timelength |
---|
| 526 | if (wprim(ilon,ilat,ilev,itim).eq.miss_val) flagfft=.false. |
---|
| 527 | enddo |
---|
| 528 | |
---|
| 529 | if (flagfft) then |
---|
| 530 | |
---|
| 531 | ! 1/ windowing to improve spectral analysis |
---|
| 532 | var(:)=wprim(ilon,ilat,ilev,:)*wndow(:) |
---|
| 533 | ! 2/ FFT computation |
---|
| 534 | !---FFTW routines |
---|
| 535 | call dfftw_execute_dft_r2c(planf,var,fftvar) |
---|
| 536 | !--- |
---|
| 537 | ! 3/ Amplitude of the FFT, for spectral analysis |
---|
| 538 | fftaw(ilon,ilat,ilev,1)=abs(fftvar(1))/M_fft |
---|
| 539 | do itim=2,M_fft |
---|
| 540 | fftaw(ilon,ilat,ilev,itim) = abs(fftvar(itim))/M_fft |
---|
| 541 | enddo |
---|
| 542 | fftaw(ilon,ilat,ilev,M_fft+1)=abs(fftvar(M_fft+1))/M_fft |
---|
| 543 | do itim=M_fft+2,timelength |
---|
| 544 | fftaw(ilon,ilat,ilev,itim) = 0. |
---|
| 545 | enddo |
---|
| 546 | |
---|
| 547 | ! 4/ filtering the FFT in three regions |
---|
| 548 | ! filtering + normalisation (low freq) |
---|
| 549 | fltfft(:) = fftvar(:)*filtrelf(:)/timelength |
---|
| 550 | ! 5/ backward FFT for each region |
---|
| 551 | !---FFTW routines |
---|
| 552 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 553 | !--- |
---|
| 554 | ! 6/ reverse the windowing |
---|
| 555 | wlf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 556 | |
---|
| 557 | ! filtering + normalisation (band freq) |
---|
| 558 | fltfft(:) = fftvar(:)*filtrebf(:)/timelength |
---|
| 559 | !---FFTW routines |
---|
| 560 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 561 | !--- |
---|
| 562 | wbf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 563 | |
---|
| 564 | ! filtering + normalisation (high freq) |
---|
| 565 | fltfft(:) = fftvar(:)*filtrehf(:)/timelength |
---|
| 566 | !---FFTW routines |
---|
| 567 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 568 | !--- |
---|
| 569 | whf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 570 | |
---|
| 571 | else |
---|
| 572 | fftaw(ilon,ilat,ilev,itim) = miss_val |
---|
| 573 | wlf(ilon,ilat,ilev,itim) = miss_val |
---|
| 574 | wbf(ilon,ilat,ilev,itim) = miss_val |
---|
| 575 | whf(ilon,ilat,ilev,itim) = miss_val |
---|
| 576 | endif ! flagfft |
---|
| 577 | |
---|
| 578 | endif !ok_out(3) |
---|
| 579 | |
---|
| 580 | ! For temperature field |
---|
| 581 | if (ok_out(4)) then |
---|
| 582 | |
---|
| 583 | flagfft=.true. |
---|
| 584 | do itim=1,timelength |
---|
| 585 | if (Tprim(ilon,ilat,ilev,itim).eq.miss_val) flagfft=.false. |
---|
| 586 | enddo |
---|
| 587 | |
---|
| 588 | if (flagfft) then |
---|
| 589 | |
---|
| 590 | ! 1/ windowing to improve spectral analysis |
---|
| 591 | var(:)=Tprim(ilon,ilat,ilev,:)*wndow(:) |
---|
| 592 | ! 2/ FFT computation |
---|
| 593 | !---FFTW routines |
---|
| 594 | call dfftw_execute_dft_r2c(planf,var,fftvar) |
---|
| 595 | !--- |
---|
| 596 | ! 3/ Amplitude of the FFT, for spectral analysis |
---|
| 597 | fftaT(ilon,ilat,ilev,1)=abs(fftvar(1))/M_fft |
---|
| 598 | do itim=2,M_fft |
---|
| 599 | fftaT(ilon,ilat,ilev,itim) = abs(fftvar(itim))/M_fft |
---|
| 600 | enddo |
---|
| 601 | fftaT(ilon,ilat,ilev,M_fft+1)=abs(fftvar(M_fft+1))/M_fft |
---|
| 602 | do itim=M_fft+2,timelength |
---|
| 603 | fftaT(ilon,ilat,ilev,itim) = 0. |
---|
| 604 | enddo |
---|
| 605 | |
---|
| 606 | ! 4/ filtering the FFT in three regions |
---|
| 607 | ! filtering + normalisation (low freq) |
---|
| 608 | fltfft(:) = fftvar(:)*filtrelf(:)/timelength |
---|
| 609 | ! 5/ backward FFT for each region |
---|
| 610 | !---FFTW routines |
---|
| 611 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 612 | !--- |
---|
| 613 | ! 6/ reverse the windowing |
---|
| 614 | Tlf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 615 | |
---|
| 616 | ! filtering + normalisation (band freq) |
---|
| 617 | fltfft(:) = fftvar(:)*filtrebf(:)/timelength |
---|
| 618 | !---FFTW routines |
---|
| 619 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 620 | !--- |
---|
| 621 | Tbf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 622 | |
---|
| 623 | ! filtering + normalisation (high freq) |
---|
| 624 | fltfft(:) = fftvar(:)*filtrehf(:)/timelength |
---|
| 625 | !---FFTW routines |
---|
| 626 | call dfftw_execute_dft_c2r(planb,fltfft,fltvar) |
---|
| 627 | !--- |
---|
| 628 | Thf(ilon,ilat,ilev,:) = fltvar(:)/wndow(:) |
---|
| 629 | |
---|
| 630 | else |
---|
| 631 | fftaT(ilon,ilat,ilev,itim) = miss_val |
---|
| 632 | Tlf(ilon,ilat,ilev,itim) = miss_val |
---|
| 633 | Tbf(ilon,ilat,ilev,itim) = miss_val |
---|
| 634 | Thf(ilon,ilat,ilev,itim) = miss_val |
---|
| 635 | endif ! flagfft |
---|
| 636 | |
---|
| 637 | endif !ok_out(4) |
---|
| 638 | |
---|
| 639 | enddo |
---|
| 640 | enddo |
---|
| 641 | enddo ! lonlength |
---|
| 642 | |
---|
| 643 | !---FFTW routines |
---|
| 644 | call dfftw_destroy_plan(planf) |
---|
| 645 | call dfftw_destroy_plan(planb) |
---|
| 646 | !--- |
---|
| 647 | |
---|
| 648 | print*,"End of computations" |
---|
| 649 | |
---|
| 650 | !=============================================================================== |
---|
| 651 | ! 3. Create output files |
---|
| 652 | !=============================================================================== |
---|
| 653 | |
---|
| 654 | ! Create output files |
---|
| 655 | if (ok_out(1)) then |
---|
| 656 | ierr=NF_CREATE(outfile1,NF_CLOBBER,outfid1) |
---|
| 657 | if (ierr.ne.NF_NOERR) then |
---|
| 658 | write(*,*)"Error: could not create file ",outfile1 |
---|
| 659 | stop |
---|
| 660 | endif |
---|
| 661 | endif !ok_out(1) |
---|
| 662 | |
---|
| 663 | if (ok_out(2)) then |
---|
| 664 | ierr=NF_CREATE(outfile2,NF_CLOBBER,outfid2) |
---|
| 665 | if (ierr.ne.NF_NOERR) then |
---|
| 666 | write(*,*)"Error: could not create file ",outfile2 |
---|
| 667 | stop |
---|
| 668 | endif |
---|
| 669 | endif !ok_out(2) |
---|
| 670 | |
---|
| 671 | if (ok_out(3)) then |
---|
| 672 | ierr=NF_CREATE(outfile3,NF_CLOBBER,outfid3) |
---|
| 673 | if (ierr.ne.NF_NOERR) then |
---|
| 674 | write(*,*)"Error: could not create file ",outfile3 |
---|
| 675 | stop |
---|
| 676 | endif |
---|
| 677 | endif !ok_out(3) |
---|
| 678 | |
---|
| 679 | if (ok_out(4)) then |
---|
| 680 | ierr=NF_CREATE(outfile4,NF_CLOBBER,outfid4) |
---|
| 681 | if (ierr.ne.NF_NOERR) then |
---|
| 682 | write(*,*)"Error: could not create file ",outfile4 |
---|
| 683 | stop |
---|
| 684 | endif |
---|
| 685 | endif !ok_out(4) |
---|
| 686 | |
---|
| 687 | !=============================================================================== |
---|
| 688 | ! 3.1. Define and write dimensions |
---|
| 689 | !=============================================================================== |
---|
| 690 | |
---|
| 691 | if (ok_out(1)) & |
---|
| 692 | call write_dim(outfid1,lonlength,latlength,altlength,timelength, & |
---|
| 693 | lon,lat,plev,time,lon_dimid1,lat_dimid1,alt_dimid1,time_dimid1) |
---|
| 694 | if (ok_out(2)) & |
---|
| 695 | call write_dim(outfid2,lonlength,latlength,altlength,timelength, & |
---|
| 696 | lon,lat,plev,time,lon_dimid2,lat_dimid2,alt_dimid2,time_dimid2) |
---|
| 697 | if (ok_out(3)) & |
---|
| 698 | call write_dim(outfid3,lonlength,latlength,altlength,timelength, & |
---|
| 699 | lon,lat,plev,time,lon_dimid3,lat_dimid3,alt_dimid3,time_dimid3) |
---|
| 700 | if (ok_out(4)) & |
---|
| 701 | call write_dim(outfid4,lonlength,latlength,altlength,timelength, & |
---|
| 702 | lon,lat,plev,time,lon_dimid4,lat_dimid4,alt_dimid4,time_dimid4) |
---|
| 703 | |
---|
| 704 | !=============================================================================== |
---|
| 705 | ! 3.2. Define and write variables |
---|
| 706 | !=============================================================================== |
---|
| 707 | |
---|
| 708 | if (ok_out(1)) then |
---|
| 709 | |
---|
| 710 | datashape4d(1)=lon_dimid1 |
---|
| 711 | datashape4d(2)=lat_dimid1 |
---|
| 712 | datashape4d(3)=alt_dimid1 |
---|
| 713 | datashape4d(4)=time_dimid1 |
---|
| 714 | datashape1d =time_dimid1 |
---|
| 715 | |
---|
| 716 | call write_var1d(outfid1,datashape1d,timelength,& |
---|
| 717 | "freq ", "FFT frequencies ","s-1 ",miss_val,& |
---|
| 718 | freq ) |
---|
| 719 | |
---|
| 720 | call write_var4d(outfid1,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 721 | "fftau ", "FFT ampl of vitu ","m s-1 ",miss_val,& |
---|
| 722 | fftau ) |
---|
| 723 | |
---|
| 724 | call write_var4d(outfid1,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 725 | "ulf ", "low freq part vitu ","m s-1 ",miss_val,& |
---|
| 726 | ulf ) |
---|
| 727 | |
---|
| 728 | call write_var4d(outfid1,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 729 | "ubf ", "band freq part vitu ","m s-1 ",miss_val,& |
---|
| 730 | ubf ) |
---|
| 731 | |
---|
| 732 | call write_var4d(outfid1,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 733 | "uhf ", "high freq part vitu ","m s-1 ",miss_val,& |
---|
| 734 | uhf ) |
---|
| 735 | endif !ok_out(1) |
---|
| 736 | |
---|
| 737 | if (ok_out(2)) then |
---|
| 738 | |
---|
| 739 | datashape4d(1)=lon_dimid2 |
---|
| 740 | datashape4d(2)=lat_dimid2 |
---|
| 741 | datashape4d(3)=alt_dimid2 |
---|
| 742 | datashape4d(4)=time_dimid2 |
---|
| 743 | datashape1d =time_dimid2 |
---|
| 744 | |
---|
| 745 | call write_var4d(outfid2,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 746 | "fftav ", "FFT ampl of vitv ","m s-1 ",miss_val,& |
---|
| 747 | fftav ) |
---|
| 748 | |
---|
| 749 | call write_var4d(outfid2,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 750 | "vlf ", "low freq part vitv ","m s-1 ",miss_val,& |
---|
| 751 | vlf ) |
---|
| 752 | |
---|
| 753 | call write_var4d(outfid2,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 754 | "vbf ", "band freq part vitv ","m s-1 ",miss_val,& |
---|
| 755 | vbf ) |
---|
| 756 | |
---|
| 757 | call write_var4d(outfid2,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 758 | "vhf ", "high freq part vitv ","m s-1 ",miss_val,& |
---|
| 759 | vhf ) |
---|
| 760 | endif !ok_out(2) |
---|
| 761 | |
---|
| 762 | if (ok_out(3)) then |
---|
| 763 | |
---|
| 764 | datashape4d(1)=lon_dimid3 |
---|
| 765 | datashape4d(2)=lat_dimid3 |
---|
| 766 | datashape4d(3)=alt_dimid3 |
---|
| 767 | datashape4d(4)=time_dimid3 |
---|
| 768 | datashape1d =time_dimid3 |
---|
| 769 | |
---|
| 770 | call write_var4d(outfid3,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 771 | "fftaw ", "FFT ampl of vitw ","Pa s-1 ",miss_val,& |
---|
| 772 | fftaw ) |
---|
| 773 | |
---|
| 774 | call write_var4d(outfid3,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 775 | "wlf ", "low freq part vitw ","Pa s-1 ",miss_val,& |
---|
| 776 | wlf ) |
---|
| 777 | |
---|
| 778 | call write_var4d(outfid3,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 779 | "wbf ", "band freq part vitw ","Pa s-1 ",miss_val,& |
---|
| 780 | wbf ) |
---|
| 781 | |
---|
| 782 | call write_var4d(outfid3,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 783 | "whf ", "high freq part vitw ","Pa s-1 ",miss_val,& |
---|
| 784 | whf ) |
---|
| 785 | endif !ok_out(3) |
---|
| 786 | |
---|
| 787 | if (ok_out(4)) then |
---|
| 788 | |
---|
| 789 | datashape4d(1)=lon_dimid4 |
---|
| 790 | datashape4d(2)=lat_dimid4 |
---|
| 791 | datashape4d(3)=alt_dimid4 |
---|
| 792 | datashape4d(4)=time_dimid4 |
---|
| 793 | datashape1d =time_dimid4 |
---|
| 794 | |
---|
| 795 | call write_var4d(outfid4,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 796 | "fftaT ", "FFT ampl of temp ","K ",miss_val,& |
---|
| 797 | fftaT ) |
---|
| 798 | |
---|
| 799 | call write_var4d(outfid4,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 800 | "tlf ", "low freq part temp ","K ",miss_val,& |
---|
| 801 | Tlf ) |
---|
| 802 | |
---|
| 803 | call write_var4d(outfid4,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 804 | "tbf ", "band freq part temp ","K ",miss_val,& |
---|
| 805 | Tbf ) |
---|
| 806 | |
---|
| 807 | call write_var4d(outfid4,datashape4d,lonlength,latlength,altlength,timelength,& |
---|
| 808 | "thf ", "high freq part temp ","K ",miss_val,& |
---|
| 809 | Thf ) |
---|
| 810 | endif !ok_out(4) |
---|
| 811 | |
---|
| 812 | !!!! Close output files |
---|
| 813 | if (ok_out(1)) then |
---|
| 814 | ierr=NF_CLOSE(outfid1) |
---|
| 815 | if (ierr.ne.NF_NOERR) write(*,*) 'Error, failed to close output file ',outfile1 |
---|
| 816 | endif !ok_out(1) |
---|
| 817 | |
---|
| 818 | if (ok_out(2)) then |
---|
| 819 | ierr=NF_CLOSE(outfid2) |
---|
| 820 | if (ierr.ne.NF_NOERR) write(*,*) 'Error, failed to close output file ',outfile2 |
---|
| 821 | endif !ok_out(2) |
---|
| 822 | |
---|
| 823 | if (ok_out(3)) then |
---|
| 824 | ierr=NF_CLOSE(outfid3) |
---|
| 825 | if (ierr.ne.NF_NOERR) write(*,*) 'Error, failed to close output file ',outfile3 |
---|
| 826 | endif !ok_out(3) |
---|
| 827 | |
---|
| 828 | if (ok_out(4)) then |
---|
| 829 | ierr=NF_CLOSE(outfid4) |
---|
| 830 | if (ierr.ne.NF_NOERR) write(*,*) 'Error, failed to close output file ',outfile4 |
---|
| 831 | endif !ok_out(4) |
---|
| 832 | |
---|
| 833 | |
---|
| 834 | end program |
---|