source: BOL/Class_Reg/geo2reg.F90 @ 5447

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

Petit oubli dans le changement du type des nbins
LF

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