da_get_var_3d_real_cdf.inc
References to this file elsewhere.
1 subroutine da_get_var_3d_real_cdf(file, var, data, &
2 i1, i2, i3, time, debug)
3
4 !-----------------------------------------------------------------------
5 ! Purpose: TBD
6 !-----------------------------------------------------------------------
7
8 implicit none
9
10 #include "netcdf.inc"
11
12 integer, intent(in) :: i1, i2, i3, time
13 character (len=80), intent(in) :: file
14 logical, intent(in) :: debug
15 character (len=*), intent(in) :: var
16 real, dimension(i1,i2,i3), intent(out) :: data
17 real(kind=8), dimension(i1,i2,i3) :: tmp
18 real(kind=4), dimension(i1,i2,i3) :: tmp4
19
20 character (len=80) :: varnam
21
22 integer :: cdfid, rcode, id_data
23 integer :: ndims, natts, idims(10), istart(10),iend(10), dimids(10)
24 integer :: i, ivtype
25
26 cdfid = ncopn(file, NCNOWRIT, rcode)
27
28 if (rcode /= 0) then
29 write(unit=stdout,fmt=*) ' error openiing netcdf file ', trim(file)
30 stop
31 end if
32
33 id_data = ncvid(cdfid, var, rcode)
34
35 rcode = nf_inq_var(cdfid, id_data, varnam, ivtype, ndims, dimids, natts)
36
37 if (debug) then
38 write(unit=*, fmt='(3a,i6)') ' get_var_3d_real_cdf: dims for ',var,' ',ndims
39 write(unit=*, fmt='(a,i6)') ' ivtype=', ivtype
40 write(unit=*, fmt='(a, a)') ' varnam=', trim(varnam)
41 write(unit=*, fmt='(a,i6)') ' kind(data)=', kind(data)
42 end if
43
44 do i=1,ndims
45 rcode = nf_inq_dimlen(cdfid, dimids(i), idims(i))
46 if (debug) write(unit=*, fmt='(a,2i6)') ' dimension ',i,idims(i)
47 end do
48
49 ! check the dimensions
50
51 if ((i1 /= idims(1)) .or. &
52 (i2 /= idims(2)) .or. &
53 (i3 /= idims(3)) .or. &
54 (time > idims(4)) ) then
55
56 write(unit=stdout,fmt=*) ' error in 3d_var_real read, dimension problem '
57 write(unit=stdout,fmt=*) i1, idims(1)
58 write(unit=stdout,fmt=*) i2, idims(2)
59 write(unit=stdout,fmt=*) i3, idims(3)
60 write(unit=stdout,fmt=*) time, idims(4)
61 write(unit=stdout,fmt=*) ' error stop '
62 stop
63 end if
64
65 ! get the data
66
67 istart(1) = 1
68 iend(1) = i1
69 istart(2) = 1
70 iend(2) = i2
71 istart(3) = 1
72 iend(3) = i3
73 istart(4) = time
74 iend(4) = 1
75
76 if ((ivtype == NF_real) .and. (kind(data) == 4)) then
77 call ncvgt(cdfid,id_data,istart,iend,data,rcode)
78 else if ((ivtype == NF_DOUBLE) .and. (kind(data) == 4)) then
79 call ncvgt(cdfid,id_data,istart,iend,tmp,rcode)
80 data = tmp
81 else if ((ivtype == NF_DOUBLE) .and. (kind(data) == 8)) then
82 call ncvgt(cdfid,id_data,istart,iend,data,rcode)
83 else if ((ivtype == NF_REAL) .and. (kind(data) == 8)) then
84 call ncvgt(cdfid,id_data,istart,iend,tmp4,rcode)
85 data = tmp4
86 else
87 write(unit=*, fmt='(a, i6)') &
88 'Unrecognizable ivtype:', ivtype
89 stop
90 end if
91
92 if (debug) then
93 write(unit=*, fmt='(a,e24.12)') ' Sample data=', data(1,1,1)
94 end if
95
96 call ncclos(cdfid,rcode)
97
98 end subroutine da_get_var_3d_real_cdf
99
100