source: trunk/mesoscale/LMD_MM_MARS/SRC/ARWpost/src/output_module.F90 @ 16

Last change on this file since 16 was 11, checked in by aslmd, 14 years ago

spiga@svn-planeto:ajoute le modele meso-echelle martien

File size: 18.8 KB
Line 
1MODULE output_module
2   
3   USE input_module
4   USE module_model_basics
5   USE module_arrays
6   
7   integer                             :: time, rec, ivars   !!! process time ; rec in .dat file ; variables in .dat file
8   integer                             :: cunit, dunit       !!! .ctl and .dat output files
9   character (len=128)                 :: ctlfile, datfile, v5dfile
10   character (len=200), dimension(200) :: ctl_var_string     !!! List of fields foe .ctl file
11   character (len=2000)                :: could_not_find     !!! Keep list of what we wrote out
12   character (len=19)                  :: tdef_date          !!! Output start date
13   character (len=40)                  :: tdef               !!! tdef string in .ctl file
14   character (len=10), dimension(100)  :: varnames
15   integer                             :: numvars
16
17   CONTAINS
18
19   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
20   ! Name: create_ctl     
21   ! Purpose: Create the .ctl file                 
22   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
23   SUBROUTINE create_ctl ()
24
25      IMPLICIT NONE
26 
27      ! Local variables
28      real, dimension(4)               :: xlat_a, xlon_a
29      real                             :: dxll, slack, temp
30      real                             :: abslatmin, abslonmax, abslonmin, abslatmax
31      integer                          :: i, k, ipoints, jpoints, ilon
32      integer                          :: ipole
33      real                             :: xi, r, re, xj, dx_ps, earthr, radpd
34      real                             :: rlat, rlong, ref_lon, ref_lat
35      integer                          :: file_start, file_end
36      integer                          :: we_dim, sn_dim
37
38
39     
40      file_start = INDEX(datfile,"/",BACK=.TRUE.)
41      file_end   = len_trim(datfile)
42      write(cunit,'("dset ^",A)') datfile(file_start+1:file_end)
43#ifdef bytesw
44   write (cunit, '("options byteswapped")')
45#endif
46      write(cunit, '("undef 1.e37")')
47
48      IF ( output_title /= '   ') THEN
49        write(cunit,'("title ",A)') trim(output_title)
50      ELSE
51        write(cunit,'("title ",A)') trim(title)
52      END IF
53
54      IF (map_proj == 1) THEN                                     !! Lambert Projection
55
56
57!!!!!!!
58!!!!!!! TEMP
59!!!!!!!
60IF ( mercator_defs ) THEN
61   WRITE(cunit,'("xdef ",i4," levels ")') west_east_dim
62     DO i = 1,west_east_dim
63        WRITE(cunit,*) XLONG(i,int(west_east_dim/2))
64     ENDDO
65   WRITE(cunit,'("ydef ",i4," levels ")') south_north_dim
66     DO i = 1,south_north_dim
67        WRITE(cunit,*) XLAT(int(south_north_dim/2),i)
68     END DO
69ELSE
70
71        ! Make sure truelat1 is the larger number
72         IF (truelat1 < truelat2) THEN
73            temp = truelat1
74            truelat1 = truelat2
75            truelat2 = temp
76         END IF
77
78         !!WRITE(cunit,'("pdef ",i4," ",i3," lcc ",f9.5," ",f10.5," ", &   !! rounding off
79         WRITE(cunit,'("pdef ",i4," ",i3," lcc ",f7.3," ",f8.3," ", &   
80&             f8.3," ",f8.3," ",f9.5," ",f9.5," ",f10.5," ",         &
81&             f10.3," ",f10.3)')                                       &
82              west_east_dim,south_north_dim,cen_lat,cen_lon,         &
83              (west_east_dim+1)/2.,(south_north_dim+1)/2.,           &
84              truelat1,truelat2,stand_lon,dx,dy
85
86
87         xlon_a(1) = XLONG(1,1)
88         xlon_a(2) = XLONG(1,south_north_dim)
89         xlon_a(3) = XLONG(west_east_dim,1)
90         xlon_a(4) = XLONG(west_east_dim,south_north_dim)
91
92         !check for dateline
93         ilon = 0
94         if ( abs(xlon_a(1) - xlon_a(3)) .GT. 180. ) ilon = 1
95         if ( abs(xlon_a(2) - xlon_a(4)) .GT. 180. ) ilon = 1
96
97         abslatmin = minval(XLAT)
98         abslatmax = maxval(XLAT)
99         abslonmin =  99999.
100         abslonmax = -99999.
101         DO i=1,4
102           IF ( xlon_a(i) .lt. 0.0 .AND. ilon .eq. 1 ) THEN
103             abslonmin=min(abslonmin,360.+xlon_a(i))
104             abslonmax=max(abslonmax,360.+xlon_a(i))
105           ELSE
106             abslonmin=min(abslonmin,xlon_a(i))
107             abslonmax=max(abslonmax,xlon_a(i))
108           ENDIF
109         END DO
110
111!         dxll = (dx/1000.)/111./2.
112         dxll = (dx/1000.)/59./2.
113         slack = ((dx/1000.)*3.)/100.
114         ipoints = int(((abslonmax-abslonmin)+(3*slack))/dxll)
115         jpoints = int(((abslatmax-abslatmin)+(3*slack))/dxll)
116
117         WRITE(cunit,'("xdef ",i4," linear ",f10.5," ",f12.8)') ipoints, &
118                      (abslonmin-1.5*slack),dxll
119         WRITE(cunit,'("ydef ",i4," linear ",f10.5," ",f12.8)') jpoints, &
120                      (abslatmin-1.5*slack),dxll
121
122!!!!!!!!!!!!!
123!!!!!!!!!!!!! TEMP
124!!!!!!!!!!!!!
125ENDIF
126
127
128      ENDIF
129
130      IF (map_proj == 2) THEN                                     !! Polar Stereo Projection
131
132         ipole = 1
133         IF (truelat1 .lt. 0.0) ipole = -1
134         radpd  = .01745329
135!         earthr = 6.3712E6
136         earthr = 3.3972E6
137
138         dx_ps = dx/( (1.0+sin(ipole*TRUELAT1*radpd))/(1.0+sin(60*radpd)) ) !! true dx (meter) at 60degrees
139         re    = ( earthr * (1.0 + sin(60*radpd)) ) / dx_ps          !!
140
141
142         IF ( ipole == 1 ) THEN
143           rlong  = (XLONG(1,1) - stand_lon) * radpd
144           rlat   = XLAT(1,1) * radpd
145
146           r    =  (re * cos(rlat)) / (1.0 + sin(rlat))
147           xi   =  (1. - r * sin(rlong))
148           xj   =  (1. + r * cos(rlong))
149           
150           write(cunit,'("pdef ",i3," ",i3," nps ",f8.3," ",f8.3," ", &
151&                f12.4," ",f12.7)') west_east_dim, south_north_dim,   &
152                 xi, xj, stand_lon, dx_ps*0.001
153
154         ELSE IF ( ipole == -1 ) THEN
155           rlong  = (180.0 + XLONG(1,1) - stand_lon) * radpd
156           rlat   = ipole*XLAT(1,1) * radpd
157
158           r    =  (re * cos(rlat)) / (1.0 + sin(rlat))
159           xi   =  (1. - ipole * r * sin(rlong))
160           xj   =  (1. + r * cos(rlong))
161
162           write(cunit,'("pdef ",i3," ",i3," sps ",f8.3," ",f8.3," ", &
163&                f12.4," ",f12.7)') west_east_dim, south_north_dim,   &
164                 xi, xj, (180+stand_lon), -0.001*dx_ps
165         END IF
166
167
168         abslonmin = minval(XLONG)
169         abslonmax = maxval(XLONG)
170         abslatmin = minval(XLAT)
171         abslatmax = maxval(XLAT)
172   
173!!         dxll = (dx_ps/1000.)/111./2.
174         dxll = (dx_ps/1000.)/59./2.
175         ipoints = int(((abslonmax-abslonmin)+2.)/dxll)-2
176         jpoints = int(((abslatmax-abslatmin)+2.)/dxll)-3
177
178         ref_lon = abslonmin-1.
179         IF ( abslonmax .ge. 175. .AND. abslonmax .le. 180.) THEN
180           IF ( abslonmin .ge. -180. .AND. abslonmin .le. -175.) THEN
181             ref_lon = 0.0
182           END IF
183         END IF
184
185         ref_lat = abslatmin-1.
186         !!IF ( ipole == -1 ) ref_lat = abslatmax+1.
187         IF ( ipole == -1 ) ref_lat = abslatmin
188   
189         WRITE(cunit,'("xdef ",i4," linear ",f6.1," ",f12.8)') ipoints, &
190                      ref_lon,dxll
191         WRITE(cunit,'("ydef ",i4," linear ",f6.1," ",f12.8)') jpoints, &
192                      ref_lat,dxll
193
194      ENDIF
195
196      IF (map_proj == 3) THEN                         !! Mercator Projection
197   
198          !check for dateline
199          ilon = 0
200          IF ( abs(XLONG(1,1) - XLONG(west_east_dim,1)) .GT. 180. ) ilon = 1
201          IF ( abs(XLONG(1,south_north_dim) - XLONG(west_east_dim,south_north_dim)) .GT. 180. ) ilon = 1
202   
203          IF ( mercator_defs .AND. map_proj == 3 ) THEN
204            WRITE(cunit,'("xdef ",i4," levels ")') west_east_dim
205            DO i = 1,west_east_dim
206              WRITE(cunit,*) XLONG(i,1)
207              !!IF ( XLONG(i,1) .lt. 0.0 ) THEN
208                !!WRITE(cunit,*) XLONG(i,1) + 180.0
209              !!ELSE
210                !!WRITE(cunit,*) XLONG(i,1)
211              !!END IF
212            END DO
213          ELSE
214            IF ( ilon == 1 ) THEN
215              WRITE(cunit,'("xdef ",i4," linear ",f9.4," ",f8.4)')   &
216                    west_east_dim,XLONG(1,1), &
217                    (abs(XLONG(1,1)-(360.+XLONG(west_east_dim,south_north_dim)))/(west_east_dim-1))
218            ELSE
219              WRITE(cunit,'("xdef ",i4," linear ",f9.4," ",f8.4)')   &
220                    west_east_dim,XLONG(1,1), &
221                    (abs(XLONG(1,1)-XLONG(west_east_dim,south_north_dim))/(west_east_dim-1))
222            ENDIF
223          ENDIF
224
225          IF ( mercator_defs .AND. map_proj == 3 ) THEN
226             WRITE(cunit,'("ydef ",i4," levels ")') south_north_dim
227             DO i = 1,south_north_dim
228               WRITE(cunit,*) XLAT(1,i)
229             END DO
230          ELSE
231             WRITE(cunit,'("ydef ",i4," linear ",f9.4," ",f8.4)')   &
232                       south_north_dim,XLAT(1,1),(abs(XLAT(1,1)-XLAT(west_east_dim,south_north_dim))/(south_north_dim-1))
233          ENDIF
234
235      ENDIF
236
237
238      IF (map_proj == 0) THEN                         !! Ideal
239        we_dim = west_east_dim
240        sn_dim = south_north_dim
241        IF (we_dim == 2) we_dim = 1
242        IF (sn_dim == 2) sn_dim = 1
243        write (cunit, '("xdef  ",i4, " linear 0 0.0001")') we_dim
244        write (cunit, '("ydef  ",i4, " linear 0 0.0001")') sn_dim
245      ENDIF
246      !! END of projection specific section
247
248      IF (iprogram .ge. 6 .and. (vertical_type=='p'.or.vertical_type=='z')) THEN
249        WRITE(cunit,'("zdef  ",i3, " levels  ")') number_of_zlevs
250        DO k = 1,number_of_zlevs
251           WRITE(cunit,'(f10.5)') interp_levels(k)
252        END DO
253      ELSE
254        WRITE(cunit,'("zdef ",i4," linear 1 1  ")') bottom_top_dim
255      ENDIF
256
257       IF (tdef_date == '0000-00-00_00:00:00') THEN
258         tdef = '   1 linear 00z01jan2000  1hr'
259       ELSE
260         IF ( output_type == 'grads' ) THEN
261                CALL tdef_calc_grads ()
262         ELSE
263                CALL tdef_calc ()
264         ENDIF
265       END IF
266       WRITE(cunit,'("tdef ",A)') tdef
267
268       WRITE(cunit,'("VARS  ",i3)') ivars
269       DO i = 1,ivars
270          WRITE(cunit,'(A)') trim(ctl_var_string(i))
271       ENDDO
272       WRITE(cunit,'("ENDVARS")')
273       
274       DO i=1,iatts
275         WRITE(cunit,'(A)') trim(catts(i))
276       ENDDO
277
278
279   END SUBROUTINE create_ctl
280
281   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
282   ! Name: tdef_calc     
283   ! Purpose: Make the tdef line for the .ctl file
284   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
285   SUBROUTINE tdef_calc ()
286
287     IMPLICIT NONE
288
289     character (len=3)              :: cmonth
290     integer                        :: year,month,day,hour
291 
292     tdef = '                               '
293 
294     !! Time comes in as YYYY-MM-DD_HH:MM:SS
295     !!                  1234 67 90 23 45 89
296     READ(tdef_date(1:4),*)   year
297     READ(tdef_date(6:7),*)   month
298     READ(tdef_date(9:10),*)  day
299     READ(tdef_date(12:13),*) hour
300
301    IF (month == 1) cmonth = 'M1'
302    IF (month == 2) cmonth = 'M2'
303    IF (month == 3) cmonth = 'M3'
304    IF (month == 4) cmonth = 'M4'
305    IF (month == 5) cmonth = 'M5'
306    IF (month == 6) cmonth = 'M6'
307    IF (month == 7) cmonth = 'M7'
308    IF (month == 8) cmonth = 'M8'
309    IF (month == 9) cmonth = 'M9'
310    IF (month ==10) cmonth = 'M10'
311    IF (month ==11) cmonth = 'M11'
312    IF (month ==12) cmonth = 'M12'
313
314    IF ( year < 100 ) year = 2000
315
316    WRITE (tdef,'(i4," linear ",i2.2,"Z",i2.2,A,i4,"  ",i6,"MN")') time, hour, day, cmonth, year, interval_seconds
317
318   END SUBROUTINE tdef_calc
319
320   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
321   ! Name: tdef_calc_grads     
322   ! Purpose: Make the tdef line for the .ctl file
323   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
324   SUBROUTINE tdef_calc_grads ()
325
326     IMPLICIT NONE
327
328     character (len=3)              :: cmonth
329     integer                        :: year,month,day,hour
330
331     tdef = '                               '
332
333     !! Time comes in as YYYY-MM-DD_HH:MM:SS
334     !!                  1234 67 90 23 45 89
335     READ(tdef_date(1:4),*)   year
336     READ(tdef_date(6:7),*)   month
337     READ(tdef_date(9:10),*)  day
338     READ(tdef_date(12:13),*) hour
339
340    IF (month == 1) cmonth = 'jan'
341    IF (month == 2) cmonth = 'fev'
342    IF (month == 3) cmonth = 'mar'
343    IF (month == 4) cmonth = 'avr'
344    IF (month == 5) cmonth = 'may'
345    IF (month == 6) cmonth = 'jun'
346    IF (month == 7) cmonth = 'jul'
347    IF (month == 8) cmonth = 'aug'
348    IF (month == 9) cmonth = 'sep'
349    IF (month ==10) cmonth = 'oct'
350    IF (month ==11) cmonth = 'nov'
351    IF (month ==12) cmonth = 'dec'
352
353    IF ( year < 100 ) year = 2000
354
355    WRITE (tdef,'(i4," linear ",i2.2,"Z",i2.2,A,i4,"  ",i6,"MN")') time, hour, day, cmonth, year, interval_seconds
356
357   END SUBROUTINE tdef_calc_grads
358
359
360   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
361   ! Name: write_dat
362   ! Purpose: Write the data to the .dat file. Also keep record of what has been
363   !          written for creation of the .ctl file
364   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
365   SUBROUTINE write_dat (data_out, nx, ny, nz, cname, cdesc, cunits)
366
367      IMPLICIT NONE
368
369      ! ------ some vis5d stuff, contain MAXVARS,MAXTIMES,MAXROWS,MAXCOLUMNS,MAXLEVELS
370      include 'v5df.h'
371
372      ! Arguments
373      integer                        :: nx, ny, nz
374      real, dimension(nx,ny,nz)      :: data_out
375      character (len=128)            :: cname, cunits, cdesc
376
377      ! Local variables
378      integer                        :: ii, jj, kk
379      integer                        :: is_there, ret, ivar
380      character (len=20)             :: dummy
381      real, dimension(ny,nx,nz)      :: v5d_data_out
382
383
384!     IF ( output_type == 'ncdf' ) THEN   !!! Create the netcdf file
385!        IF ( nx == 2 ) nx = 1
386!        IF ( ny == 2 ) ny = 1
387!        !DO kk=1,nz
388!        !  rec = rec + 1
389!        !  WRITE(dunit,rec=rec) ((data_out(ii,jj,kk),ii=1,nx),jj=1,ny)
390!        !END DO
391!        CALL writencdf(unite nx,ny,nz,"u","W.m-2",3,data_out)
392!     END IF
393
394
395
396      IF ( output_type == 'grads' ) THEN   !!! Create the grads file
397        IF ( nx == 2 ) nx = 1
398        IF ( ny == 2 ) ny = 1
399        DO kk=1,nz
400          rec = rec + 1
401          WRITE(dunit,rec=rec) ((data_out(ii,jj,kk),ii=1,nx),jj=1,ny)
402        END DO
403      END IF
404
405!!!!****MARS IDL
406      IF ( output_type == 'idl' ) THEN   !!! Create the IDL ascii file
407        IF ( nx == 2 ) nx = 1
408        IF ( ny == 2 ) ny = 1
409        DO kk=1,nz
410          rec = rec + 1
411!!!****MARS
412!          WRITE(dunit,rec=rec) ((data_out(ii,jj,kk),ii=1,nx),jj=1,ny)
413          WRITE(dunit,FMT='(F20.8)') ((data_out(ii,jj,kk),ii=1,nx),jj=1,ny)
414        END DO
415      END IF
416
417
418
419#ifdef V5D
420      IF ( output_type == 'v5d' ) THEN   !!! Create the vis5d file
421        DO kk = 1,nz
422        DO jj = 1,nx
423        DO ii = 1,ny
424          v5d_data_out(ii,jj,kk) = data_out(jj,ny-ii+1,kk)
425        END DO
426        END DO
427        END DO
428
429
430        ivar = 0
431        DO ii = 1,numvars
432          IF ( trim(varnames(ii)) == trim(cname) ) THEN
433            ivar = ii
434          END IF
435        END DO
436        ret = v5dwrite ( time+1, ivar, v5d_data_out )
437      END IF
438#endif
439
440
441      IF (time == 0) THEN
442        IF ( cname(1:3) == 'clf' ) then
443          dummy = ",clfr,"
444        ELSE
445          dummy = ","//trim(cname)//","
446        END IF
447        is_there = INDEX(could_not_find,trim(dummy))
448        DO WHILE ( is_there /= 0 )
449          could_not_find = could_not_find(1:is_there)//could_not_find(is_there+len_trim(dummy):len_trim(could_not_find))
450          is_there = INDEX(could_not_find,trim(dummy))
451        END DO
452        ivars = ivars + 1
453        WRITE(ctl_var_string(ivars),'(A,"  ",i4,"  0  ",A," (",A,")")') cname(1:10),nz,trim(cdesc),trim(cunits)
454      ENDIF
455
456
457   END SUBROUTINE write_dat
458
459!   SUBROUTINE writencdf(unite,nx,ny,nz,"u","W.m-2",3,data_out)
460!
461!
462!        ! NETCDF file
463!        fichnom="diagfi.nc"
464!        if (firstnom.eq.'1234567890') then
465!        firstnom = nom
466!        ierr = NF_CREATE  (fichnom, NF_CLOBBER, nid)
467!        ! Build dimensions
468!        ierr = NF_DEF_DIM (nid, "Time",      NF_UNLIMITED, idim)
469!        ierr = NF_DEF_DIM (nid, "latitude",  ny,           idim_lat)
470!        ierr = NF_DEF_DIM (nid, "longitude", nx,           idim_lon)
471!        ierr = NF_DEF_DIM (nid, "altitude",  nz,           idim_alt)
472!        ! Write dimensions attributes
473!        ierr = NF_DEF_VAR (nid, "Time", NF_FLOAT, 1, idim, varid)
474!        ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name",  4, "Time")
475!        ierr = NF_PUT_ATT_TEXT (nid, varid, 'units',     29, "Time")
476!        ierr = NF_ENDDEF(nid)
477!        !
478!        ierr = NF_DEF_VAR (nid, "latitude",  NF_FLOAT,  1, idim_lat, varid)
479!        ierr = NF_PUT_ATT_TEXT (nid, varid, 'units',     13, "degrees_north")
480!        ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name", 14, "North latitude")
481!        ierr = NF_ENDDEF(nid)
482!        ierr = NF_DEF_VAR (nid, "longitude", NF_FLOAT,  1, idim_lon, varid)
483!        ierr = NF_PUT_ATT_TEXT (nid, varid, 'units',     12, "degrees_east")
484!        ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name", 14, "East longitude")
485!        ierr = NF_ENDDEF(nid)
486!        ierr = NF_DEF_VAR (nid, "altitude",  NF_FLOAT,  1, idim_alt, varid)
487!        ierr = NF_PUT_ATT_TEXT (nid, varid, 'units',      2, "km")
488!        ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name", 10, "pseudo-alt")
489!        ierr = NF_PUT_ATT_TEXT (nid, varid, 'positive',   2, "up")
490!        ierr = NF_ENDDEF(nid)
491!        else
492!        ierr = NF_OPEN(fichnom,NF_WRITE,nid)
493!        endif
494!        ! Advance time
495!        ierr = NF_INQ_VARID (nid, "Time", varid)
496!        ierr = NF_PUT_VARA_REAL (nid, varid, ntime, 1, date)
497!        ! Write field
498!        ierr = NF_INQ_VARID(nid,nom,varid)
499!        if (ierr /= NF_NOERR) then
500!        ! Create variable if does not exist
501!        ierr = NF_INQ_DIMID(nid,"longitude",id(1))
502!        ierr = NF_INQ_DIMID(nid,"latitude",id(2))
503!        ierr = NF_INQ_DIMID(nid,"altitude",id(3))
504!        ierr = NF_INQ_DIMID(nid,"Time",id(4))
505!        write (*,*) "====================================="
506!        write (*,*) "NETCDF: create ",nom
507!        call def_var (nid,nom,titre,unite,4,id,varid,ierr)
508!        endif
509!   END SUBROUTINE writencdf
510
511
512!!==================================================================================
513!!==================================================================================
514!
515!subroutine def_var(nid,name,title,units,nbdim,dim,nvarid,ierr)
516!
517!implicit none
518!
519!include "netcdf.inc"
520!
521!character (len=*) :: title,units,name
522!integer :: nid,nbdim,nvarid,ierr
523!integer, dimension(nbdim) :: dim
524!
525!ierr=NF_REDEF(nid)
526!#ifdef NC_DOUBLE
527!ierr = NF_DEF_VAR (nid,adjustl(name),NF_DOUBLE,nbdim,dim,nvarid)
528!#else
529!ierr = NF_DEF_VAR (nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
530!#endif
531!if(ierr/=NF_NOERR) then
532!   write(*,*) NF_STRERROR(ierr)
533!   stop "in def_var"
534!endif
535!ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title",
536!len_trim(adjustl(title)),adjustl(title))
537!if(ierr/=NF_NOERR) then
538!   write(*,*) NF_STRERROR(ierr)
539!   stop "in def_var"
540!endif
541!ierr = NF_PUT_ATT_TEXT (nid, nvarid, "units",
542!len_trim(adjustl(units)),adjustl(units))
543!if(ierr/=NF_NOERR) then
544!   write(*,*) NF_STRERROR(ierr)
545!   stop "in def_var"
546!endif
547!ierr = NF_ENDDEF(nid)
548!
549!end subroutine def_var
550
551
552END MODULE output_module
Note: See TracBrowser for help on using the repository browser.