source: codes/icosagcm/trunk/src/write_field.f90 @ 15

Last change on this file since 15 was 15, checked in by ymipsl, 12 years ago

Update on 3D dynamic

YM

File size: 18.6 KB
Line 
1module write_field
2USE genmod
3implicit none
4
5  integer, parameter :: MaxWriteField = 1000
6  integer, dimension(MaxWriteField),save :: FieldId
7  integer, dimension(MaxWriteField),save :: FieldVarId
8  integer, dimension(MaxWriteField),save :: FieldIndex
9  character(len=255), dimension(MaxWriteField) ::  FieldName 
10   
11  integer,save :: NbField = 0
12 
13  contains
14 
15    function GetFieldIndex(name)
16    implicit none
17      integer          :: GetFieldindex
18      character(len=*) :: name
19   
20      character(len=255) :: TrueName
21      integer            :: i
22       
23     
24      TrueName=TRIM(ADJUSTL(name))
25   
26      GetFieldIndex=-1
27      do i=1,NbField
28        if (TrueName==FieldName(i)) then
29          GetFieldIndex=i
30          exit
31        endif
32      enddo
33    end function GetFieldIndex
34 
35
36     
37    subroutine WriteField_gen(name,Field,dimx,dimy,dimz)
38    implicit none
39!    include 'netcdf.inc'
40      character(len=*) :: name
41      integer :: dimx,dimy,dimz
42      real,dimension(dimx,dimy,dimz) :: Field
43      integer,dimension(dimx*dimy*dimz) :: ndex
44      integer :: status
45      integer :: index
46      integer :: start(4)
47      integer :: count(4)
48     
49           
50      Index=GetFieldIndex(name)
51      if (Index==-1) then
52!        call CreateNewField(name,dimx,dimy,dimz)
53        Index=GetFieldIndex(name)
54      else
55        FieldIndex(Index)=FieldIndex(Index)+1.
56      endif
57     
58      start(1)=1
59      start(2)=1
60      start(3)=1
61      start(4)=FieldIndex(Index)
62
63      count(1)=dimx
64      count(2)=dimy
65      count(3)=dimz
66      count(4)=1
67
68!      status = NF_PUT_VARA_DOUBLE(FieldId(Index),FieldVarId(Index),start,count,Field)
69!      status = NF_SYNC(FieldId(Index))
70     
71    end subroutine WriteField_gen
72
73
74    SUBROUTINE Writefield(name_in,field,nind)
75    USE netcdf
76    USE domain_mod
77    use field_mod
78    USE dimensions
79    USE geometry
80    IMPLICIT NONE
81      CHARACTER(LEN=*),INTENT(IN) :: name_in
82      TYPE(t_field),POINTER :: field(:)
83      INTEGER,OPTIONAL,INTENT(IN) :: nind
84      REAL(r8),ALLOCATABLE :: field_val2d(:)
85      REAL(r8),ALLOCATABLE :: field_val3d(:,:)
86      REAL(r8),ALLOCATABLE :: field_val4d(:,:,:)
87      TYPE(t_domain),POINTER :: d
88      INTEGER :: Index
89      INTEGER :: ind,i,j,k,n,ncell
90      INTEGER :: iie,jje,iin,jjn
91      INTEGER :: status
92      CHARACTER(len=255) :: name
93      CHARACTER(len=255) :: str_ind
94      INTEGER :: ind_b,ind_e
95      INTEGER :: halo_size
96      LOGICAL :: single
97     
98     
99      name=TRIM(ADJUSTL(name_in))
100
101      IF (PRESENT(nind)) THEN
102        name=TRIM(name)//"_"//TRIM(int2str(nind))
103        PRINT *,"NAME",nind,int2str(nind),name
104        ind_b=nind
105        ind_e=nind
106        halo_size=1
107        single=.TRUE.
108      ELSE
109        ind_b=1
110        ind_e=ndomain
111        halo_size=0
112        single=.FALSE.
113      ENDIF     
114
115      Index=GetFieldIndex(name)
116      if (Index==-1) then
117        call create_header(name,field,nind)
118        Index=GetFieldIndex(name)
119      else
120        FieldIndex(Index)=FieldIndex(Index)+1.
121      endif
122     
123      IF (Field(ind_b)%field_type==field_T) THEN
124        ncell=1
125        DO ind=ind_b,ind_e
126          d=>domain(ind)
127          IF (Field(ind)%field_type/=field_T) THEN
128            PRINT *,"Writefield, grille non geree"
129            RETURN
130          ENDIF
131
132        n=0
133        DO j=d%jj_begin-halo_size,d%jj_end+halo_size
134          DO i=d%ii_begin-halo_size,d%ii_end+halo_size
135            IF (d%own(i,j) .OR. single) n=n+1
136          ENDDO
137        ENDDO
138
139        IF (field(ind)%ndim==2) THEN
140          ALLOCATE(Field_val2d(n))
141          n=0
142          DO j=d%jj_begin-halo_size,d%jj_end+halo_size
143            DO i=d%ii_begin-halo_size,d%ii_end+halo_size
144              k=d%iim*(j-1)+i
145              IF (d%own(i,j) .OR. single) THEN
146                n=n+1
147                Field_val2d(n)=field(ind)%rval2d(k)
148              ENDIF
149            ENDDO
150          ENDDO
151          status=NF90_PUT_VAR(FieldId(Index),FieldVarId(index),Field_val2d,start=(/ ncell,FieldIndex(Index) /),count=(/n,1 /))
152          DEALLOCATE(field_val2d)
153        ELSE IF (field(ind)%ndim==3) THEN
154          ALLOCATE(Field_val3d(n,size(field(ind)%rval3d,2)))
155          n=0
156          DO j=d%jj_begin-halo_size,d%jj_end+halo_size
157            DO i=d%ii_begin-halo_size,d%ii_end+halo_size
158              k=d%iim*(j-1)+i
159              IF (d%own(i,j) .OR. single) THEN
160                n=n+1
161                Field_val3d(n,:)=field(ind)%rval3d(k,:)
162              ENDIF
163            ENDDO
164          ENDDO
165           status=NF90_PUT_VAR(FieldId(Index),FieldVarId(index),Field_val3d,start=(/ ncell,1,FieldIndex(Index) /),   &
166                               count=(/n,size(field(1)%rval3d,2),1 /))
167          DEALLOCATE(field_val3d)
168        ELSE IF (field(1)%ndim==4) THEN
169          ALLOCATE(Field_val4d(n,size(field(ind)%rval4d,2),size(field(ind)%rval4d,3)))
170          n=0
171          DO j=d%jj_begin-halo_size,d%jj_end+halo_size
172            DO i=d%ii_begin-halo_size,d%ii_end+halo_size
173              k=d%iim*(j-1)+i
174              IF (d%own(i,j) .OR. single) THEN
175                n=n+1
176                Field_val4d(n,:,:)=field(ind)%rval4d(k,:,:)
177              ENDIF
178            ENDDO
179          ENDDO
180
181          status=NF90_PUT_VAR(FieldId(Index),FieldVarId(index),Field_val4d,start=(/ ncell,1,1,FieldIndex(Index) /),   &
182                              count=(/n,size(field(1)%rval4d,2),size(field(1)%rval4d,3),1 /))
183          DEALLOCATE(field_val4d)
184        ENDIF
185       
186!        PRINT *,NF90_STRERROR(status)
187        ncell=ncell+n
188
189     ENDDO
190     
191     ELSE IF (Field(ind_b)%field_type==field_Z) THEN
192        ncell=1
193        n=0
194        DO ind=ind_b,ind_e
195          d=>domain(ind)
196          CALL swap_geometry(ind)
197          CALL swap_dimensions(ind)
198 
199          n=0
200          DO j=jj_begin+1,jj_end
201            DO i=ii_begin,ii_end-1
202              n=n+1
203            ENDDO
204          ENDDO
205
206          DO j=jj_begin,jj_end-1
207            DO i=ii_begin+1,ii_end
208              n=n+1
209            ENDDO
210          ENDDO
211
212        IF (field(ind)%ndim==2) THEN
213          ALLOCATE(Field_val2d(n))
214
215          n=0
216          DO j=jj_begin+1,jj_end
217            DO i=ii_begin,ii_end-1
218              n=n+1
219              k=iim*(j-1)+i
220              Field_val2d(n)=field(ind)%rval2d(k+z_down)
221            ENDDO
222          ENDDO
223
224          DO j=jj_begin,jj_end-1
225            DO i=ii_begin+1,ii_end
226              n=n+1
227              k=iim*(j-1)+i
228              Field_val2d(n)=field(ind)%rval2d(k+z_up)
229            ENDDO
230          ENDDO
231
232          status=NF90_PUT_VAR(FieldId(Index),FieldVarId(index),Field_val2d,start=(/ ncell,FieldIndex(Index) /),count=(/n,1 /))
233          DEALLOCATE(field_val2d)
234
235        ELSE IF (field(ind)%ndim==3) THEN
236          ALLOCATE(Field_val3d(n,size(field(ind)%rval3d,2)))
237          n=0
238          DO j=jj_begin+1,jj_end
239            DO i=ii_begin,ii_end-1
240              n=n+1
241              k=iim*(j-1)+i
242              Field_val3d(n,:)=field(ind)%rval3d(k+z_down,:)
243            ENDDO
244          ENDDO
245
246          DO j=jj_begin,jj_end-1
247            DO i=ii_begin+1,ii_end
248              n=n+1
249              k=iim*(j-1)+i
250              Field_val3d(n,:)=field(ind)%rval3d(k+z_up,:)
251            ENDDO
252          ENDDO
253           status=NF90_PUT_VAR(FieldId(Index),FieldVarId(index),Field_val3d,start=(/ ncell,1,FieldIndex(Index) /),   &
254                               count=(/n,size(field(1)%rval3d,2),1 /))
255          DEALLOCATE(field_val3d)
256        ELSE IF (field(1)%ndim==4) THEN
257          ALLOCATE(Field_val4d(n,size(field(ind)%rval4d,2),size(field(ind)%rval4d,3)))
258          n=0
259          DO j=jj_begin+1,jj_end
260            DO i=ii_begin,ii_end-1
261              n=n+1
262              k=iim*(j-1)+i
263              Field_val4d(n,:,:)=field(ind)%rval4d(k+z_down,:,:)
264            ENDDO
265          ENDDO
266
267          DO j=jj_begin,jj_end-1
268            DO i=ii_begin+1,ii_end
269              n=n+1
270              k=iim*(j-1)+i
271              Field_val4d(n,:,:)=field(ind)%rval4d(k+z_up,:,:)
272            ENDDO
273          ENDDO
274
275          status=NF90_PUT_VAR(FieldId(Index),FieldVarId(index),Field_val4d,start=(/ ncell,1,1,FieldIndex(Index) /),   &
276                              count=(/n,size(field(1)%rval4d,2),size(field(1)%rval4d,3),1 /))
277          DEALLOCATE(field_val4d)
278        ENDIF
279       
280!        PRINT *,NF90_STRERROR(status)
281        ncell=ncell+n
282
283     ENDDO
284     
285     ENDIF
286     status=NF90_SYNC(FieldId(Index))
287     
288    END SUBROUTINE Writefield
289     
290           
291    SUBROUTINE Create_header(name,field,nind)
292    USE netcdf
293    USE field_mod
294    USE domain_mod
295    USE spherical_geom_mod
296    USE dimensions
297    USE geometry
298    IMPLICIT NONE
299      CHARACTER(LEN=*) :: name
300      TYPE(t_field),POINTER :: field(:)
301      INTEGER,OPTIONAL,INTENT(IN) :: nind
302      INTEGER :: ncell
303      INTEGER :: nvert
304      REAL(rstd),ALLOCATABLE :: lon(:),lat(:),bounds_lon(:,:),bounds_lat(:,:)
305      TYPE(t_domain),POINTER :: d
306      INTEGER :: nvertId,ncid,lonId,latId,bounds_lonId,bounds_latId,timeId,ncellId
307      INTEGER :: dim3id,dim4id
308      INTEGER :: status
309      INTEGER :: ind,i,j,k,n
310      INTEGER :: iie,jje,iin,jjn
311      INTEGER :: ind_b,ind_e
312      INTEGER :: halo_size
313      LOGICAL :: single 
314      INTEGER :: nij
315         
316      NbField=NbField+1
317      FieldName(NbField)=TRIM(ADJUSTL(name))
318      FieldIndex(NbField)=1
319     
320      IF (PRESENT(nind)) THEN
321        ind_b=nind
322        ind_e=nind
323        halo_size=1
324        single=.TRUE.
325      ELSE
326        ind_b=1
327        ind_e=ndomain
328        halo_size=0
329        single=.FALSE.
330      ENDIF
331     
332      ncell=0
333     
334      IF (Field(ind_b)%field_type==field_T) THEN
335        nvert=6
336       
337        DO ind=ind_b,ind_e
338          d=>domain(ind)
339       
340          DO j=d%jj_begin-halo_size,d%jj_end+halo_size
341            DO i=d%ii_begin-halo_size,d%ii_end+halo_size
342              IF (single .OR. domain(ind)%own(i,j)) ncell=ncell+1
343            ENDDO
344          ENDDO
345
346        END DO
347     
348        status = NF90_CREATE(TRIM(ADJUSTL(name))//'.nc', NF90_CLOBBER, ncid)
349        FieldId(NbField)=ncid
350        status = NF90_DEF_DIM(ncid,'cell',ncell,ncellId)
351        status = NF90_DEF_DIM(ncid,'nvert',nvert,nvertid)
352
353        IF (Field(ind_b)%ndim==3)  THEN
354          status = NF90_DEF_DIM(ncid,'Z',size(field(ind_b)%rval3d,2),dim3id)
355        ELSE IF (Field(1)%ndim==4) THEN
356          status = NF90_DEF_DIM(ncid,'Z',size(field(ind_b)%rval4d,2),dim3id)
357          status = NF90_DEF_DIM(ncid,'Q',size(field(ind_b)%rval4d,3),dim4id)
358        ENDIF
359     
360        status = NF90_DEF_DIM(ncid,'time_counter',NF90_UNLIMITED,timeId)
361     
362        status = NF90_DEF_VAR(ncid,'lon',NF90_DOUBLE,(/ ncellId /),lonId)
363        status = NF90_PUT_ATT(ncid,lonId,"long_name","longitude")
364        status = NF90_PUT_ATT(ncid,lonId,"units","degrees_east")
365        status = NF90_PUT_ATT(ncid,lonId,"bounds","bounds_lon")
366        status = NF90_DEF_VAR(ncid,'lat',NF90_DOUBLE,(/ ncellId /),latId)
367        status = NF90_PUT_ATT(ncid,latId,"long_name","latitude")
368        status = NF90_PUT_ATT(ncid,latId,"units","degrees_north")
369        status = NF90_PUT_ATT(ncid,latId,"bounds","bounds_lat")
370        status = NF90_DEF_VAR(ncid,'bounds_lon',NF90_DOUBLE,(/ nvertid,ncellId /),bounds_lonId)
371        status = NF90_DEF_VAR(ncid,'bounds_lat',NF90_DOUBLE,(/ nvertid,ncellId /),bounds_latId)
372
373        IF (Field(ind_b)%ndim==2) THEN
374          status = NF90_DEF_VAR(ncid,TRIM(ADJUSTL(name)),NF90_DOUBLE,(/ ncellId,timeId /),FieldVarId(NbField))
375        ELSE IF (Field(ind_b)%ndim==3) THEN
376          status = NF90_DEF_VAR(ncid,TRIM(ADJUSTL(name)),NF90_DOUBLE,(/ ncellId,dim3id,timeId /),FieldVarId(NbField))
377        ELSE IF (Field(ind_b)%ndim==4) THEN
378          status = NF90_DEF_VAR(ncid,TRIM(ADJUSTL(name)),NF90_DOUBLE,(/ ncellId,dim3id,dim4id,timeId /),FieldVarId(NbField))
379        ENDIF
380 
381        status = NF90_PUT_ATT(ncid,FieldVarId(NbField),"coordinates","lon lat")
382       
383        status = NF90_ENDDEF(ncid)     
384
385        ncell=1
386        DO ind=ind_b,ind_e
387          d=>domain(ind)
388 
389          n=0
390          DO j=d%jj_begin-halo_size,d%jj_end+halo_size
391            DO i=d%ii_begin-halo_size,d%ii_end+halo_size
392              IF (single .OR. d%own(i,j)) n=n+1
393            ENDDO
394          ENDDO
395         
396         ALLOCATE(lon(n),lat(n),bounds_lon(0:nvert-1,n),bounds_lat(0:nvert-1,n))
397         
398          n=0 
399          DO j=d%jj_begin-halo_size,d%jj_end+halo_size
400            DO i=d%ii_begin-halo_size,d%ii_end+halo_size
401                IF (d%own(i,j) .OR. single) THEN
402                n=n+1
403                CALL xyz2lonlat(d%xyz(:,i,j),lon(n),lat(n))
404                lon(n)=lon(n)*180/Pi
405                lat(n)=lat(n)*180/Pi
406                DO k=0,5
407                  CALL xyz2lonlat(d%vertex(:,k,i,j),bounds_lon(k,n), bounds_lat(k,n))
408                  bounds_lat(k,n)=bounds_lat(k,n)*180/Pi
409                  bounds_lon(k,n)=bounds_lon(k,n)*180/Pi
410                ENDDO
411              ENDIF
412            ENDDO
413          ENDDO
414          status=NF90_PUT_VAR(ncid,lonid,REAL(lon,r8),start=(/ ncell /),count=(/ n /))
415          status=NF90_PUT_VAR(ncid,latid,REAL(lat,r8),start=(/ ncell /),count=(/ n /))
416          status=NF90_PUT_VAR(ncid,bounds_lonId,REAL(bounds_lon,r8),start=(/ 1,ncell /),count=(/ nvert,n /))
417          status=NF90_PUT_VAR(ncid,bounds_latId,REAL(bounds_lat,r8),start=(/ 1,ncell /),count=(/ nvert,n /))
418 
419          ncell=ncell+n
420          DEALLOCATE(lon,lat,bounds_lon,bounds_lat)
421      END DO
422
423    ELSE IF (Field(ind_b)%field_type==field_Z) THEN
424        nvert=3
425        DO ind=ind_b,ind_e
426          d=>domain(ind)
427       
428          DO j=d%jj_begin+1,d%jj_end
429            DO i=d%ii_begin,d%ii_end-1
430              ncell=ncell+1
431            ENDDO
432          ENDDO
433
434          DO j=d%jj_begin,d%jj_end-1
435            DO i=d%ii_begin+1,d%ii_end
436              ncell=ncell+1
437            ENDDO
438          ENDDO
439
440        END DO
441     
442        status = NF90_CREATE(TRIM(ADJUSTL(name))//'.nc', NF90_CLOBBER, ncid)
443        FieldId(NbField)=ncid
444        status = NF90_DEF_DIM(ncid,'cell',ncell,ncellId)
445        status = NF90_DEF_DIM(ncid,'nvert',nvert,nvertid)
446
447        IF (Field(ind_b)%ndim==3)  THEN
448          status = NF90_DEF_DIM(ncid,'Z',size(field(ind_b)%rval3d,2),dim3id)
449        ELSE IF (Field(1)%ndim==4) THEN
450          status = NF90_DEF_DIM(ncid,'Z',size(field(ind_b)%rval4d,2),dim3id)
451          status = NF90_DEF_DIM(ncid,'Q',size(field(ind_b)%rval4d,3),dim4id)
452        ENDIF
453     
454        status = NF90_DEF_DIM(ncid,'time_counter',NF90_UNLIMITED,timeId)
455     
456        status = NF90_DEF_VAR(ncid,'lon',NF90_DOUBLE,(/ ncellId /),lonId)
457        status = NF90_PUT_ATT(ncid,lonId,"long_name","longitude")
458        status = NF90_PUT_ATT(ncid,lonId,"units","degrees_east")
459        status = NF90_PUT_ATT(ncid,lonId,"bounds","bounds_lon")
460        status = NF90_DEF_VAR(ncid,'lat',NF90_DOUBLE,(/ ncellId /),latId)
461        status = NF90_PUT_ATT(ncid,latId,"long_name","latitude")
462        status = NF90_PUT_ATT(ncid,latId,"units","degrees_north")
463        status = NF90_PUT_ATT(ncid,latId,"bounds","bounds_lat")
464        status = NF90_DEF_VAR(ncid,'bounds_lon',NF90_DOUBLE,(/ nvertid,ncellId /),bounds_lonId)
465        status = NF90_DEF_VAR(ncid,'bounds_lat',NF90_DOUBLE,(/ nvertid,ncellId /),bounds_latId)
466
467        IF (Field(ind_b)%ndim==2) THEN
468          status = NF90_DEF_VAR(ncid,TRIM(ADJUSTL(name)),NF90_DOUBLE,(/ ncellId,timeId /),FieldVarId(NbField))
469        ELSE IF (Field(ind_b)%ndim==3) THEN
470          status = NF90_DEF_VAR(ncid,TRIM(ADJUSTL(name)),NF90_DOUBLE,(/ ncellId,dim3id,timeId /),FieldVarId(NbField))
471        ELSE IF (Field(ind_b)%ndim==4) THEN
472          status = NF90_DEF_VAR(ncid,TRIM(ADJUSTL(name)),NF90_DOUBLE,(/ ncellId,dim3id,dim4id,timeId /),FieldVarId(NbField))
473        ENDIF
474 
475        status = NF90_PUT_ATT(ncid,FieldVarId(NbField),"coordinates","lon lat")
476       
477        status = NF90_ENDDEF(ncid)     
478
479        ncell=1
480        DO ind=ind_b,ind_e
481          d=>domain(ind)
482          CALL swap_geometry(ind)
483          CALL swap_dimensions(ind)
484 
485          n=0
486          DO j=jj_begin+1,jj_end
487            DO i=ii_begin,ii_end-1
488              n=n+1
489            ENDDO
490          ENDDO
491
492          DO j=jj_begin,jj_end-1
493            DO i=ii_begin+1,ii_end
494              n=n+1
495            ENDDO
496          ENDDO
497
498         ALLOCATE(lon(n),lat(n),bounds_lon(0:nvert-1,n),bounds_lat(0:nvert-1,n))
499         
500          n=0 
501       
502          DO j=jj_begin+1,jj_end
503            DO i=ii_begin,ii_end-1
504              nij=(j-1)*iim+i
505              n=n+1
506              CALL xyz2lonlat(xyz_v(nij+z_down,:)/radius,lon(n),lat(n))
507              lon(n)=lon(n)*180/Pi
508 !             IF (lon(n)<0) lon(n)=lon(n)+360
509              lat(n)=lat(n)*180/Pi
510              CALL xyz2lonlat(xyz_i(nij,:)/radius,bounds_lon(0,n), bounds_lat(0,n))
511              CALL xyz2lonlat(xyz_i(nij+t_ldown,:)/radius,bounds_lon(1,n), bounds_lat(1,n))
512              CALL xyz2lonlat(xyz_i(nij+t_rdown,:)/radius,bounds_lon(2,n), bounds_lat(2,n))
513
514              DO k=0,2
515                bounds_lat(k,n)=bounds_lat(k,n)*180/Pi
516                bounds_lon(k,n)=bounds_lon(k,n)*180/Pi
517!                IF (bounds_lon(k,n)<0) bounds_lon(k,n)=bounds_lon(k,n)+360
518              ENDDO
519            ENDDO
520          ENDDO
521
522          DO j=jj_begin,jj_end-1
523            DO i=ii_begin+1,ii_end
524              nij=(j-1)*iim+i
525              n=n+1
526              CALL xyz2lonlat(xyz_v(nij+z_up,:)/radius,lon(n),lat(n))
527              lon(n)=lon(n)*180/Pi
528!              IF (lon(n)<0) lon(n)=lon(n)+360
529              lat(n)=lat(n)*180/Pi
530              CALL xyz2lonlat(xyz_i(nij,:)/radius,bounds_lon(0,n), bounds_lat(0,n))
531              CALL xyz2lonlat(xyz_i(nij+t_rup,:)/radius,bounds_lon(1,n), bounds_lat(1,n))
532              CALL xyz2lonlat(xyz_i(nij+t_lup,:)/radius,bounds_lon(2,n), bounds_lat(2,n))
533
534              DO k=0,2
535                bounds_lat(k,n)=bounds_lat(k,n)*180/Pi
536                bounds_lon(k,n)=bounds_lon(k,n)*180/Pi
537!                IF (bounds_lon(k,n)<0) bounds_lon(k,n)=bounds_lon(k,n)+360
538              ENDDO
539            ENDDO
540          ENDDO
541         
542         
543          status=NF90_PUT_VAR(ncid,lonid,REAL(lon,r8),start=(/ ncell /),count=(/ n /))
544          status=NF90_PUT_VAR(ncid,latid,REAL(lat,r8),start=(/ ncell /),count=(/ n /))
545          status=NF90_PUT_VAR(ncid,bounds_lonId,REAL(bounds_lon,r8),start=(/ 1,ncell /),count=(/ nvert,n /))
546          status=NF90_PUT_VAR(ncid,bounds_latId,REAL(bounds_lat,r8),start=(/ 1,ncell /),count=(/ nvert,n /))
547 
548          ncell=ncell+n
549          DEALLOCATE(lon,lat,bounds_lon,bounds_lat)
550      END DO         
551    ENDIF
552
553
554     
555!      status = NF90_CLOSE(ncid)
556
557    end subroutine Create_Header
558   
559 
560  function int2str(int)
561    implicit none
562    integer, parameter :: MaxLen=10
563    integer,intent(in) :: int
564    character(len=MaxLen) :: int2str
565    logical :: flag
566    integer :: i
567    flag=.true.
568   
569    i=int
570   
571    int2str=''
572    do while (flag)
573      int2str=CHAR(MOD(i,10)+48)//int2str
574      i=i/10
575      if (i==0) flag=.false.
576    enddo
577  end function int2str
578
579end module write_field
580 
Note: See TracBrowser for help on using the repository browser.