indexing typo in phys/module_fr_sfire_atm.F/interpolate_wind2fire_height/interpolate_h
[wrf-fire.git] / other / netcdf_crop / netcdf_crop.F
blob024d24c631c3fca947d429cbab9b005baaacf8b3
2 program netcdf_crop
3   use netcdf
4   implicit none
6   integer n,i,j
7   character (len=256) :: carg,fin,fout
8   character (len=nf90_max_name) nname,tname,aname
9   integer tval,ncidin,ncidout,ndims,nvars,ngattr,nvattr,nunlim,nvdims
10   integer ilen,itmp,nlen,tid
11   integer,dimension(nf90_max_var_dims) :: dimids,odimids
12   integer,dimension(nf90_max_dims) :: dlens,dconv
13   integer,dimension(nf90_max_vars) :: varids,xtype
14   integer,dimension(:),allocatable :: ivals
15   character,dimension(:),allocatable :: cvals
16   real(kind=kind(0.0d0)),dimension(:),allocatable :: fvals
17   integer,dimension(:),allocatable :: start,count,ostart,ocount
19   n=iargc()
20   if(n.ne.4)then
21     print*,'usage: netcdf_crop infilename outfilename timename time'
22     print*,'infilename: netcdf file to import'
23     print*,'outfilename: netcdf file to export'
24     print*,'timename: name of time dimension in netcdf file'
25     print*,'time: time to export'
26     goto 999
27   endif
28   call getarg(1,fin)
29   call getarg(2,fout)
30   call getarg(4,carg)
31   call getarg(3,tname)
32   read(carg,*)tval
34   call check(nf90_open(fin,nf90_nowrite,ncidin))
35   call check(nf90_create(fout,nf90_clobber,ncidout))
36   call check(nf90_inquire(ncidin,nDimensions=ndims, &
37                                  nVariables=nvars, &
38                                  nAttributes=ngattr, &
39                                  unlimitedDimId=nunlim))
40   call check(nf90_inq_dimid(ncidin,tname,tid))
41   do i=1,ndims
42     call check(nf90_inquire_dimension(ncidin,i,name=nname,len=dlens(i)))
43     if(i.ne.tid)then
44       call check(nf90_def_dim(ncidout,nname,dlens(i),dconv(i)))
45     else
46       call check(nf90_def_dim(ncidout,nname,nf90_unlimited,dconv(i)))
47     endif
48   enddo
49   do i=1,ngattr
50     call check(nf90_inq_attname(ncidin,nf90_global,i,nname))
51 !    call check(nf90_inquire_attribute(ncidin,nf90_global,nname))
52     call check(nf90_copy_att(ncidin,nf90_global,nname,ncidout,nf90_global))
53   enddo
54   do i=1,nvars
55     call check(nf90_inquire_variable(ncidin,i,name=nname,&
56                                               xtype=xtype(i),&
57                                               ndims=nvdims,&
58                                               dimids=dimids,&
59                                               nAtts=nvattr))
60     call convert_dimids(dconv,nvdims,dimids,odimids)
61     call check(nf90_def_var(ncidout,nname,xtype(i),odimids(1:nvdims),varids(i)))
62     do j=1,nvattr
63       call check(nf90_inq_attname(ncidin,i,j,aname))
64       call check(nf90_copy_att(ncidin,i,aname,ncidout,varids(i)))
65     enddo
66   enddo
67   call check(nf90_enddef(ncidout))
68   do i=1,nvars
69     call check(nf90_inquire_variable(ncidin,i,dimids=dimids,ndims=nvdims))
70     allocate(start(1:nvdims),count(1:nvdims),ostart(1:nvdims),ocount(1:nvdims))
71     do j=1,nvdims
72       if(dimids(j).ne.tid)then
73         start(j)=1
74         count(j)=dlens(dimids(j))
75         ostart(j)=1
76         ocount(j)=dlens(dimids(j))
77       else
78         start(j)=tval
79         count(j)=1
80         ostart(j)=1
81         ocount(j)=1
82       endif
83     enddo
84     call get_nlen(dlens,nvdims,dimids,tid,nlen)
85     if(xtype(i).eq.nf90_char)then
86       allocate(cvals(1:nlen))
87       call check(nf90_get_var(ncidin,i,cvals,start=start,count=count))
88       call check(nf90_put_var(ncidout,varids(i),cvals,start=ostart,&
89                                                       count=ocount))
90       deallocate(cvals)
91     elseif(xtype(i).eq.nf90_byte.or.&
92            xtype(i).eq.nf90_short.or.&
93            xtype(i).eq.nf90_int)then
94       allocate(ivals(1:nlen))
95       call check(nf90_get_var(ncidin,i,ivals,start=start,count=count))
96       call check(nf90_put_var(ncidout,varids(i),ivals,start=ostart,&
97                                                       count=ocount))
98       deallocate(ivals)
99     elseif(xtype(i).eq.nf90_float.or.&
100            xtype(i).eq.nf90_double)then
101       allocate(fvals(1:nlen))
102       call check(nf90_get_var(ncidin,i,fvals,start=start,count=count))
103       call check(nf90_put_var(ncidout,varids(i),fvals,start=ostart,&
104                                                       count=ocount))
105       deallocate(fvals)
106     else
107       call check(nf90_ebadtype)
108     endif
109     deallocate(start,count,ostart,ocount)
110   enddo
112   call check(nf90_close(ncidin))
113   call check(nf90_close(ncidout))
114   999 continue
115 contains
116   subroutine check(status)
117     integer, intent ( in) :: status
118     
119     if(status /= nf90_noerr) then 
120       print *, trim(nf90_strerror(status))
121       stop "Stopped"
122     end if
123   end subroutine check  
125   subroutine convert_dimids(dconv,nd,din,dout)
126   use netcdf
127   implicit none
128   integer,dimension(nf90_max_dims),intent(in) :: dconv
129   integer,intent(in)::nd
130   integer,dimension(nd),intent(in)::din
131   integer,dimension(nd),intent(out)::dout
132   integer i
134   do i=1,nd
135     dout(i)=dconv(din(i))
136   enddo 
138   end subroutine convert_dimids
140   subroutine get_nlen(dlens,nd,dimids,tid,nlen)
141   use netcdf
142   implicit none
143   integer,dimension(nf90_max_dims),intent(in)::dlens
144   integer,intent(in)::nd
145   integer,dimension(nd),intent(in)::dimids
146   integer,intent(in)::tid
147   integer,intent(out)::nlen
148   integer::i
149   nlen=1
150   do i=1,nd
151     if(dimids(i).ne.tid)then
152       nlen=nlen*dlens(dimids(i))
153     endif
154   enddo
155   end subroutine get_nlen
157 end program netcdf_crop