source: BOL/trunk/Class_Reg/geo2reg.F90 @ 611

Last change on this file since 611 was 611, checked in by Laurent Fairhead, 20 years ago

Pb avec la lecture du fichier config et un deallocate qui manquait
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 18.0 KB
Line 
1!
2! $Header$
3!
4!-------------------------------------------------------------------
5! Make relationship between any variable and w500 at resolution 2.5x2.5
6! -> mean seasonal cycle of the relationship over the period 1985-1989
7! -> interannual anomalies (by removing the mean seasonal cycle)
8! -> decomposition of interannual anomalies into thermo, dynam and mix components
9!
10! Specify:
11!       - undefined value of the data
12!       - min/max/step of the w500 bins
13!       - nb of months in timeseries
14!       - montly index of January 1985
15!
16! Sandrine Bony, March 2004
17! Repris pour faire du netcdf, L. Fairhead 2005/02
18!-------------------------------------------------------------------
19
20  implicit none
21
22#include "netcdf.inc"
23
24  integer imax,jmax,an,mmax,nbinmax,noc
25  parameter (imax=360,jmax=180,mmax=1700,an=12,nbinmax=200,noc=5)
26
27  integer nrun
28  parameter (nrun=1)
29
30  real zmsk(imax,jmax)
31  real surfmax,zsurf(imax,jmax)
32  logical ok_msk(imax,jmax)
33
34  character*20 namegcm, run
35
36  real surftest
37  real lat0, undef_x
38!  parameter (lat0 = 30.0) ! limits tropical belt
39  parameter (undef_x = 1.e+20 ) ! undefined value
40
41  real, allocatable, dimension(:,:) :: surf, msk, w, x
42  real, allocatable, dimension(:,:) :: x_binw
43  real, allocatable, dimension(:,:) :: w_binw
44  real, allocatable, dimension(:,:) :: xpdf,nx_binw,sx_binw
45  real, allocatable, dimension(:)   :: xpdfmean
46
47  integer nb_bin
48  real min_bin, max_bin, step_bin, w_mult, x_mult
49  real xpdftot(mmax)
50  real xpdfmeantot
51
52
53
54  INTEGER i,j,m,ir,im,jm,lm,itime,n, il, l
55  real undef, undef_w, pi, x1, w1, msk1
56  integer :: iostat
57  character (len=512)        :: line_read
58  character (len=20)         :: first_part
59  character (len=128)       :: second_part 
60  character cnat*4, bintype*4
61
62! netcdf stuff
63  real         :: time
64  integer, dimension(:), allocatable :: var_dim
65  INTEGER      :: ierr, nvdims, lunout, in_file_id, ierrs, ndims, natts
66  integer      :: lonid, latid, in_time_id, varid, out_file_id, ngatts, levid
67  integer      :: nbin_id, var_bin_id, time_id, var_time_id, in_var_time_id
68  INTEGER      :: var_pdf_id, var_nx_id, var_w_id, var_x_id, outlev_id
69  integer      :: var_lev_id
70  integer, dimension(4)  :: vdims, start, count
71  character (len=80) :: long_name, varname, attname
72  character (len=132) :: in_file
73  real, dimension(:), allocatable ::lat
74
75  logical :: var_3d = .false.
76
77  integer :: iargc
78  external iargc 
79!---------------------------------------------------------------------
80!-- General specifications :
81!---------------------------------------------------------------------
82
83!******************************* FIN INTERFACE *********************
84! Lecture du fichier de configuration
85  open (20, IOSTAT=iostat, file='config.def',form='formatted')
86  if (iostat /= 0) THEN
87    write(*,*)'Erreur ouverture du fichier config.def'
88    stop
89  endif
90
91  config: do
92    read(20,'(A)',iostat=iostat)line_read
93    if (iostat /= 0) exit
94    line_read = trim(line_read)
95    IF (INDEX(line_read, '#') /= 1) THEN
96      first_part = trim(line_read(1:INDEX(line_read, '=')-1))
97      second_part = trim(line_read(INDEX(line_read, '=')+1:))
98      selectcase(first_part)
99        case('bintype')
100          bintype = trim(second_part)
101        case('cnat')
102          cnat = trim(second_part)
103        case('lunout')
104          read(second_part,'(i)') lunout
105        case('lat0')
106          read(second_part,'(f)') lat0
107      end select
108    endif
109  enddo config
110  if (iostat > 0) then
111    write(lunout,*) &
112 &  'Probleme de lecture du fichier config.def, iostat = ',iostat
113    stop
114  endif
115  close(20)
116 
117! undefined value for the output:
118  undef  = 999.999 ! for output only
119
120!------------------------------------------------------------
121! Definition of the w500 bins:
122!------------------------------------------------------------
123
124  if (bintype.eq.'w500') then
125    min_bin = -200.0 
126    max_bin =  200.0
127    step_bin =  10.0
128    step_bin =  5.0
129    w_mult =  864.0 ! Pa/s -> hPa/day
130!   w_mult =  1.0 ! Pa/s -> hPa/day
131  endif
132
133! scale factor for surface temperature, precip and variable x:
134  x_mult =  1.0
135
136  undef_w = undef
137
138!---------------------------------------------------------------------
139!-- preliminaries:
140!---------------------------------------------------------------------
141
142  if (cnat.ne.'ocea'.and.cnat.ne.'land'.and.cnat.ne.'glob') then
143    write(*,*) 'erreur cnat ',cnat
144    stop
145  endif
146
147  nb_bin = (max_bin-min_bin)/step_bin + 1
148  if (nb_bin.gt.nbinmax) then
149    write(*,*) 'augmenter le nb de bins'
150    write(*,*) 'nbinmax, nb_bin = ',nbinmax,nb_bin
151    stop
152  endif
153
154! Il faut deux arguments a l appel: fichier et variable
155  CALL getarg(1, in_file)
156  if (iargc() == 0 .OR. iargc() /= 2 ) THEN
157    write(lunout,*)' '
158    write(lunout,*)' Utilisation de ce programme: '
159    write(lunout,*)' ./geo2reg nom_de_fichier [variable]'
160    write(lunout,*)                                    &
161    &  '        ou nom_de_fichier est le nom du fichier a traiter'
162    write(lunout,*)                                             &
163    &  '        et variable la variable a traiter [optionel]'
164       write(lunout,*)' '
165    write(lunout,*)' ./geo2reg -h sort ce message'
166    write(lunout,*)' '
167    stop
168  endif
169  CALL getarg(2, varname)
170 
171! Ouverture du fichier a traiter (histmth)
172  ierr = NF_OPEN(in_file, NF_NOwrite, in_file_id)
173  if (ierr /= NF_NOERR) then
174    write(lunout,*)NF_STRERROR(ierr)
175    stop
176  endif
177!
178! lire im, jm, lm et temps
179  ierrs = 0
180  ierr = NF_INQ_DIMID(in_file_id, 'lon', lonid) 
181  ierrs = ierrs + ierr
182  ierr = NF_INQ_DIMLEN(in_file_id, lonid, im)
183  ierrs = ierrs + ierr
184  if (ierrs /= 2 * NF_NOERR) THEN
185    write(lunout,*)'Pb. avec la lecture de la longitude'
186    stop
187  endif
188  ierrs = 0
189  ierr = NF_INQ_DIMID(in_file_id, 'lat', latid) 
190  ierrs = ierrs + ierr
191  ierr = NF_INQ_DIMLEN(in_file_id, latid, jm)
192  ierrs = ierrs + ierr
193  allocate (lat(jm))
194  ierr = NF_INQ_VARID(in_file_id,'lat', latid)
195  ierrs = ierrs + ierr
196  ierr = NF_GET_VAR_REAL(in_file_id, latid, lat)
197  ierrs = ierrs + ierr
198  if (ierrs /= 4 * NF_NOERR) THEN
199    write(lunout,*)'Pb. avec la lecture de la latitude'
200    stop
201  endif
202  ierrs = 0
203  ierr = NF_INQ_DIMID(in_file_id, 'presnivs', levid)
204  ierrs = ierrs + ierr
205  ierr = NF_INQ_DIMLEN(in_file_id, levid, lm)
206  ierrs = ierrs + ierr
207  if (ierrs /= 2 * NF_NOERR) THEN
208    write(lunout,*)'Pb. avec la lecture des niveaux verticaux'
209    stop
210  endif
211 
212
213  ierrs = 0
214  ierr = NF_INQ_DIMID(in_file_id, 'time_counter', in_time_id)
215  ierrs = ierrs + ierr
216  ierr = NF_INQ_DIMLEN(in_file_id, in_time_id, itime)
217  ierrs = ierrs + ierr
218  if (ierrs /= 2 * NF_NOERR) THEN
219    write(lunout,*)'Pb. avec la lecture du temps'
220    stop
221  endif
222
223! lecture de l aire des mailles
224  allocate (surf(im,jm))
225  ierrs = 0
226  ierr = NF_INQ_VARID(in_file_id, 'aire', varid)
227  ierrs = ierrs + ierr
228  ierr = NF_GET_VAR_REAL(in_file_id, varid, surf)
229  ierrs = ierrs + ierr
230  if (ierrs /= 2 * NF_NOERR) THEN
231    write(lunout,*)'Pb. avec la lecture de l''aire'
232    stop
233  endif
234 
235!-- compute gridbox areas:
236
237
238  pi = acos(-1.)
239  surfmax=0.
240  surftest=0.0
241  do i = 1, im
242    do j = 1, jm
243      surftest=surftest+surf(i,j)
244      if(surf(i,j).gt.surfmax) surfmax=surf(i,j)
245    enddo
246  enddo
247  write(*,*) 'sfce totale = 4pi ',surftest/pi
248
249!---------------------------------------------------------------
250! Lecture du masque
251!
252  allocate (msk(im,jm))
253  if (cnat.ne.'glob') then
254    ierrs = 0
255    ierr = NF_INQ_VARID(in_file_id, 'pourc_ter', varid)
256    ierrs = ierrs + ierr
257    ierr = NF_GET_VAR_REAL(in_file_id, varid, msk)
258    ierrs = ierrs + ierr
259    if (ierrs /= 2 * NF_NOERR) THEN
260      write(lunout,*)'Pb. avec la lecture du masque'
261      stop
262    endif
263  else
264    do i = 1, im
265      do j = 1, jm
266        msk(i,j) = undef
267      enddo
268    enddo
269  endif
270
271
272  do i = 1, im
273    do j = 1, jm
274      ok_msk(i,j) = .FALSE.
275      if (cnat.eq.'ocea' .and. msk(i,j).le.0.5) ok_msk(i,j) = .TRUE.
276      if (cnat.eq.'land' .and. msk(i,j).gt.0.5) ok_msk(i,j) = .TRUE.
277      if (cnat.eq.'glob') ok_msk(i,j) = .TRUE.
278    enddo
279  enddo
280
281!
282! Pour savoir si la variable a traiter est 3D
283!
284  ierrs = 0
285  ierr = NF_INQ_VARID(in_file_id, varname, varid)
286  ierrs = ierrs + ierr
287  ierrs = NF_INQ_VARNDIMS(in_file_id, varid, ndims)
288  if (ierrs /= 2 * NF_NOERR) THEN
289    write(lunout,*)'Pb. avec la lecture de ', varname
290    stop
291  endif
292  if (ndims == 4) var_3d = .true.   
293
294!--------------------------------------------------------------------
295! creation du fichier de sortie et entete
296!--------------------------------------------------------------------
297
298
299  ierr = NF_CREATE('out.nc', NF_CLOBBER, out_file_id)
300  if (ierr /= NF_NOERR) then
301    write(lunout,*)NF_STRERROR(ierr)
302    stop
303  endif
304  ierr = NF_INQ_NATTS(in_file_id, ngatts)
305  do i = 1, ngatts
306    ierr = NF_INQ_ATTNAME(in_file_id, NF_GLOBAL, i, attname)
307    ierr = NF_COPY_ATT(in_file_id, NF_GLOBAL, attname, out_file_id, NF_GLOBAL)
308  enddo
309!
310! definition des dimensions (nbin, level, temps)
311! nbin:
312  ierrs = 0
313  ierr = NF_DEF_DIM(out_file_id, 'nbin', nb_bin, nbin_id)
314  ierrs = ierrs + ierr
315  nvdims = 1
316  vdims(1) = nbin_id
317  ierr = NF_DEF_VAR(out_file_id, 'nbin', NF_INT, nvdims, vdims, var_bin_id)
318  ierrs = ierrs + ierr
319  ierr = NF_PUT_ATT_TEXT(out_file_id, var_bin_id, 'units', 1, '-')
320  ierrs = ierrs + ierr
321  ierr = NF_PUT_ATT_TEXT(out_file_id, var_bin_id, 'long_name', 15, &
322 &       'number of bins')
323  ierrs = ierrs + ierr
324  if (ierrs /= 4 * NF_NOERR) THEN
325    write(lunout,*)'Pb. dans la definition de nbin'
326    stop
327  endif
328!
329! presnivs
330  ierrs = 0
331  ierr = NF_DEF_DIM(out_file_id, 'presnivs', lm, outlev_id)
332  ierrs = ierrs + ierr
333  nvdims = 1
334  vdims(1) = outlev_id
335  ierr = NF_DEF_VAR(out_file_id,'presnivs',NF_FLOAT, nvdims, vdims, var_lev_id)
336  ierrs = ierrs + ierr
337  ierr = NF_INQ_VARID(in_file_id, 'presnivs', levid)
338  ierrs = ierrs + ierr
339  ierr = NF_INQ_VARNATTS(in_file_id, levid, natts)
340  ierrs = ierrs + ierr
341  do i = 1, natts
342    attname = ''
343    ierr = NF_INQ_ATTNAME(in_file_id, levid, i, attname)
344    ierrs = ierrs + ierr
345    ierr = NF_COPY_ATT(in_file_id, levid, attname, out_file_id, var_lev_id)
346    ierrs = ierrs + ierr
347  enddo
348  IF (ierrs /= (4 + natts) * NF_NOERR) THEN
349    write(lunout,*)'Pb. dans la definition de presnivs'
350    stop
351  endif
352 
353!
354! temps:
355  ierrs = 0
356  ierr = NF_DEF_DIM(out_file_id, 'time_counter', NF_UNLIMITED, time_id)
357  ierrs = ierrs + ierr
358  nvdims = 1
359  vdims(1) = time_id
360  ierr = NF_DEF_VAR(out_file_id, 'time_counter', NF_FLOAT, nvdims, vdims, &
361 &                  var_time_id)
362  ierrs = ierrs + ierr
363  ierr = NF_INQ_VARID(in_file_id, 'time_counter', in_var_time_id)
364  if (ierr /= NF_NOERR) &
365 &   ierr = NF_INQ_VARID(in_file_id, 't_ave_02592000', in_var_time_id)
366  if (ierr == NF_NOERR) then
367    ierr = NF_COPY_ATT(in_file_id, in_var_time_id, 'units', out_file_id, &
368 &                     var_time_id)
369  else
370    ierr = NF_PUT_ATT_TEXT(out_file_id, var_time_id, 'units', 34, &
371 &                         'seconds since 1860-01-01 00:00:00')
372    in_var_time_id = 0
373  endif
374  ierrs = ierrs + ierr
375  ierr = NF_PUT_ATT_TEXT(out_file_id, var_time_id, 'long_name', 9, &
376 &       'Time axis')
377  ierrs = ierrs + ierr
378  if (ierrs /= 4 * NF_NOERR) THEN
379    write(lunout,*)'Pb. dans la definition de time_counter'
380    stop
381  endif
382!
383! Definition des variables a ecrire:
384! pdf:
385  ierrs = 0
386  nvdims = 2
387  vdims(1) = nbin_id
388  vdims(2) = time_id
389  ierr = NF_DEF_VAR(out_file_id, 'pdf', NF_FLOAT, nvdims, vdims, &
390  &                 var_pdf_id)
391  ierrs = ierrs + ierr
392  ierr = NF_PUT_ATT_TEXT(out_file_id, var_pdf_id, 'long_name', 16, &
393 &       'monthly PDF w500')
394  ierrs = ierrs + ierr
395  ierr = NF_PUT_ATT_REAL(out_file_id, var_pdf_id, 'missing_value', NF_FLOAT, &
396 &       1, undef)
397  ierrs = ierrs + ierr
398  if (ierrs /= 3 * NF_NOERR) THEN
399    write(lunout,*)'Pb. dans la definition de pdf'
400    stop
401  endif
402! nx:
403  ierrs = 0
404  nvdims = 2
405  vdims(1) = nbin_id
406  vdims(2) = time_id
407  ierr = NF_DEF_VAR(out_file_id, 'nx', NF_FLOAT, nvdims, vdims, &
408  &                 var_nx_id)
409  ierrs = ierrs + ierr
410  ierr = NF_PUT_ATT_TEXT(out_file_id, var_nx_id, 'long_name', 16, &
411 &       'nb of points in w500 bin')
412  ierrs = ierrs + ierr
413  if (ierrs /= 2 * NF_NOERR) THEN
414    write(lunout,*)'Pb. dans la definition de nx'
415    stop
416  endif
417!
418! w:
419  ierrs = 0
420  nvdims = 2
421  vdims(1) = nbin_id
422  vdims(2) = time_id
423  ierr = NF_DEF_VAR(out_file_id, 'w', NF_FLOAT, nvdims, vdims, &
424  &                 var_w_id)
425  ierrs = ierrs + ierr
426  long_name = 'mean relationship w(binw)'
427  ierr = NF_PUT_ATT_TEXT(out_file_id, var_w_id, 'long_name', 25, &
428 &       long_name)
429  ierrs = ierrs + ierr
430  ierr = NF_PUT_ATT_REAL(out_file_id, var_w_id, 'missing_value', NF_FLOAT, &
431 &       1, undef)
432  ierrs = ierrs + ierr
433  if (ierrs /= 3 * NF_NOERR) THEN
434    write(lunout,*)'Pb. dans la definition de w'
435    stop
436  endif
437!
438! x:
439  ierrs = 0
440  nvdims = 2
441  vdims(1) = nbin_id
442  vdims(2) = time_id
443  if (var_3d) then
444    nvdims = 3
445    vdims(2) = outlev_id
446    vdims(3) = time_id
447  endif
448  ierr = NF_DEF_VAR(out_file_id, varname, NF_FLOAT, nvdims, vdims, &
449  &                 var_x_id)
450  ierrs = ierrs + ierr
451  long_name = 'mean relationship '//trim(varname)
452  ierr = NF_PUT_ATT_TEXT(out_file_id, var_x_id, 'long_name', len(long_name), &
453 &       long_name)
454  ierrs = ierrs + ierr
455  ierr = NF_PUT_ATT_REAL(out_file_id, var_x_id, 'missing_value', NF_FLOAT, &
456 &       1, undef)
457  ierrs = ierrs + ierr
458  if (ierrs /= 3 * NF_NOERR) THEN
459    write(lunout,*)'Pb. dans la definition de ',varname
460    stop
461  endif
462  ierr = NF_ENDDEF(out_file_id)
463
464!
465! On passe aux calculs
466!
467! initialisations:
468
469
470  allocate(xpdfmean(nb_bin))
471  allocate(xpdf(nb_bin, itime))
472  allocate(nx_binw(nb_bin, itime))
473  allocate(sx_binw(nb_bin, itime))
474  allocate(w_binw(nb_bin, itime))
475  allocate(x_binw(nb_bin, itime))
476
477   
478   
479
480
481!-- temporal loop:
482  il = 1
483  if (var_3d) then
484    il = lm
485  endif
486
487  do l = 1, il   
488  xpdfmeantot = 0.
489  xpdfmean = 0.0
490  x_binw = 0.0
491  xpdf = 0.0
492  nx_binw = 0.0
493  sx_binw = 0.0
494  w_binw = 0.0
495  xpdftot = 0.
496  do m = 1, itime ! loop over timeperiod of x_data
497
498! read data:
499    start(1) = 1
500    start(2) = 1
501    start(3) = m
502    count(1) = im
503    count(2) = jm
504    count(3) = 1
505! vit. verticale à 500hP
506    if (.not. allocated(w)) allocate (w(im,jm))
507    ierrs = 0
508    ierr = NF_INQ_VARID(in_file_id, 'w500', varid)
509    ierrs = ierrs + ierr
510    ierr = NF_GET_VARA_REAL(in_file_id, varid, start, count, w)
511    ierrs = ierrs + ierr
512    if (ierrs /= 2* NF_NOERR) THEN
513      write(lunout,*)'Pb. avec la lecture de la vit. vert. à 500hP'
514      stop
515    endif
516!
517! variable a traiter
518    if (.not. allocated(x)) allocate (x(im,jm))
519    if (var_3d) then
520      start(3) = l
521      count(3) = 1
522      start(4) = m
523      count(4) = 1
524    endif
525    ierrs = 0
526    ierr = NF_INQ_VARID(in_file_id, varname, varid)
527    ierrs = ierrs + ierr
528    ierr = NF_GET_VARA_REAL(in_file_id, varid, start, count, x)
529    ierrs = ierrs + ierr
530    if (ierrs /= 2 * NF_NOERR) THEN
531      write(lunout,*)'Pb. avec la lecture de ', varname
532      stop
533    endif
534!
535    do j=1,jm
536      do i=1,im
537        zmsk(i,j)=0.
538      enddo
539    enddo
540! tropical belt:
541    do j = 1, jm
542      if (ABS(lat(j)).le.lat0) then
543! loop over grid points:
544        do i = 1, im
545          if (x(i,j).ne.undef_x .and. w(i,j).ne.undef_w     &
546  &           .and. ok_msk(i,j) ) then
547            zmsk(i,j)=1.
548            w1 = w(i,j)*w_mult
549            x1 = x(i,j)*x_mult
550! bin w500:
551            ir = INT((w1-min_bin)/step_bin) + 1
552! monthly PDF :
553            xpdfmeantot = xpdfmeantot + surf(i,j)
554            xpdfmean(ir) = xpdfmean(ir) + surf(i,j)
555            xpdftot(m) = xpdftot(m) + surf(i,j)
556            xpdf(ir,m) = xpdf(ir,m) + surf(i,j)
557! monthly x-w relationship:
558            nx_binw(ir,m) = nx_binw(ir,m) + 1.0
559            sx_binw(ir,m) = sx_binw(ir,m) + surf(i,j)
560            x_binw(ir,m) = x_binw(ir,m) + x1*surf(i,j)
561            w_binw(ir,m) = w_binw(ir,m) + w1*surf(i,j)
562
563          endif
564        enddo ! i
565      endif ! lat
566    enddo ! j
567  enddo ! m
568
569! normalize PDF:
570
571  do ir = 1, nb_bin
572    xpdfmean(ir) = xpdfmean(ir)/xpdfmeantot
573  enddo
574
575  do m = 1, itime
576    if (xpdftot(m).gt.0.) then
577      do ir = 1, nb_bin
578        xpdf(ir,m) = xpdf(ir,m)/xpdftot(m)
579      enddo
580    else
581      do ir = 1, nb_bin
582        xpdf(ir,m) = undef
583      enddo
584    endif
585  enddo
586
587
588  do ir = 1, nb_bin
589    do m = 1, itime
590      if (nx_binw(ir,m).gt.1.) then
591!      if (nx_binw(ir,m).gt.0.) then
592        x_binw(ir,m) = x_binw(ir,m)/sx_binw(ir,m)
593        w_binw(ir,m) = w_binw(ir,m)/sx_binw(ir,m)
594      else
595        x_binw(ir,m) = undef
596        w_binw(ir,m) = undef
597      endif
598    enddo
599  enddo ! ir
600
601!
602! ecriture du fichier sortie
603!
604
605  allocate(var_dim(nb_bin))
606  do ir = 1, nb_bin
607    var_dim(ir) = min_bin + step_bin * (ir -1)
608  enddo
609  ierr = NF_PUT_VAR_INT(out_file_id, var_bin_id, var_dim)
610  if (ierr /= NF_NOERR) then
611    write(lunout,*)NF_STRERROR(ierr)
612    stop
613  endif
614  deallocate(var_dim)
615!
616! ecriture niveaux verticaux
617!
618  allocate(var_dim(lm))
619  ierrs = 0
620  ierr = NF_GET_VAR_REAL(in_file_id, levid, var_dim)
621  ierrs = ierrs + ierr
622  ierr = NF_PUT_VAR_REAL(out_file_id, var_lev_id, var_dim)
623  ierrs = ierrs + ierr
624  if (ierrs /= 2 * NF_NOERR) THEN
625    write(lunout,*)'Pb. d''ecriture niveaux verticaux'
626    stop
627  endif
628 
629 
630   
631  do m = 1, itime
632    start(1) = 1
633    start(2) = m
634    count(1) = nb_bin
635    count(2) = 1
636    ierrs = 0
637!
638! ecriture temps:
639    if (in_var_time_id == 0) then
640      time = (m - 1) * 86400. * 30.
641    else
642      ierr = NF_GET_VAR1_REAL(in_file_id, in_var_time_id, m, time)
643    endif
644    ierr = NF_PUT_VAR1_REAL(out_file_id, var_time_id, m, time)
645    ierrs = ierrs + ierr
646!
647! pdf, nb_bin, x, w:
648    ierr = NF_PUT_VARA_REAL(out_file_id, var_pdf_id, start, count, &
649   &                        xpdf(1,m))
650    ierrs = ierrs + ierr
651    ierr = NF_PUT_VARA_REAL(out_file_id, var_nx_id, start, count, &
652   &                        nx_binw(1,m))
653    ierrs = ierrs + ierr
654    ierr = NF_PUT_VARA_REAL(out_file_id, var_w_id, start, count, &
655   &                        w_binw(1,m))
656    ierrs = ierrs + ierr
657    if (var_3d) then
658      start(2) = l
659      count(2) = 1
660      start(3) = m
661      count(4) = 1
662    endif
663    ierr = NF_PUT_VARA_REAL(out_file_id, var_x_id, start, count, &
664   &                        x_binw(1,m))
665    ierrs = ierrs + ierr
666  if (ierrs /= 5 * NF_NOERR) THEN
667    write(lunout,*)'Pb. d''ecriture'
668    stop
669  endif
670  enddo ! m
671  enddo ! l
672
673  ierr = NF_CLOSE(out_file_id)
674
675  end
Note: See TracBrowser for help on using the repository browser.