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

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

Installation de la boite a outil permettant de faire du classement en regimes
a partir des fichiers histmth.nc:

geo2reg.F90 programme fortran faisant la classification
geo2reg.script un exemple de script
config.def comme son nom implique

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