Loading...
Searching...
No Matches
Functions/Subroutines
gtdata_netcdf_internal Module Reference

Functions/Subroutines

integer function, public vtable_add (var, entry)
 
integer function, public vtable_delete (var)
 
integer function, public vtable_lookup (var, entry)
 
integer function, public vtable_set_attrid (var, attrid)
 

Function/Subroutine Documentation

◆ vtable_add()

integer function, public gtdata_netcdf_internal::vtable_add ( type(gd_nc_variable), intent(out)  var,
type(gd_nc_variable_search), intent(in)  entry 
)

Definition at line 33 of file gtdata_netcdf_internal.f90.

34 type(GD_NC_VARIABLE), intent(out):: var
35 type(GD_NC_VARIABLE_SEARCH), intent(in):: entry
36 type(GD_NC_VARIABLE_ENTRY), allocatable:: tmp_table(:)
37 integer:: i, n
38
39 ! --- 必要なら初期確保 ---
40 if (.not. allocated(gdnctab)) then
41 allocate(gdnctab(gdnctab_init_size), stat=result)
42 if (result /= 0) goto 999
43 do, i = 1, gdnctab_init_size
44 gdnctab(i)%fileid = 0
45 gdnctab(i)%varid = 0
46 gdnctab(i)%dimid = 0
47 gdnctab(i)%attrid = 0
48 nullify(gdnctab(i)%dimids)
49 enddo
50 endif
51 ! --- 同じ内容が既登録ならばそれを返す (attrid は変更しない) ---
52 do, i = 1, size(gdnctab)
53 if (gdnctab(i)%fileid == entry%fileid &
54 & .and. gdnctab(i)%varid == entry%varid &
55 & .and. gdnctab(i)%dimid == entry%dimid) then
56 var = gd_nc_variable(i)
57 result = nf90_noerr
58 call dbgmessage('gtdata_netcdf_internal.add: found %d', i=(/i/))
59 return
60 endif
61 enddo
62 !
63 ! --- 空き地があればそこに割り当て ---
64 var = gd_nc_variable(-1)
65 do, i = 1, size(gdnctab)
66 if (gdnctab(i)%fileid == 0) then
67 var = gd_nc_variable(i)
68 exit
69 endif
70 enddo
71 if (var%id == -1) then
72 ! --- 空き地はなかったのだから倍幅確保 ---
73 n = size(gdnctab)
74 allocate(tmp_table(n), stat=result)
75 if (result /= 0) goto 999
76 tmp_table(:) = gdnctab(:)
77 deallocate(gdnctab, stat=result)
78 if (result /= 0) goto 999
79 allocate(gdnctab(n * 2), stat=result)
80 if (result /= 0) goto 999
81 gdnctab(1:n) = tmp_table(1:n)
82 deallocate(tmp_table, stat=result)
83 if (result /= 0) goto 999
84 !
85 gdnctab(n+2)%fileid = 0
86 gdnctab(n+2)%varid = 0
87 gdnctab(n+2)%dimid = 0
88 gdnctab(n+2)%attrid = 0
89 nullify(gdnctab(n+2)%dimids)
90 gdnctab(n+3: n*2) = gdnctab(n+2)
91 ! 確保域の先頭を使用
92 var = gd_nc_variable(n + 1)
93 endif
94 gdnctab(var%id)%fileid = entry%fileid
95 gdnctab(var%id)%varid = entry%varid
96 gdnctab(var%id)%dimid = entry%dimid
97 !
98 ! --- 次元表の確保 ---
99 call internal_build_dimids(gdnctab(var%id), result)
100 if (result /= nf90_noerr) goto 999
101 !
102 result = nf90_noerr
103 call dbgmessage('gtdata_netcdf_internal.add: added %d', i=(/var%id/))
104 return
105 !
106999 continue
107 var = gd_nc_variable(-1)
108 result = nf90_enomem
109 return
110
111 contains
112
113 subroutine internal_build_dimids(ent, stat)
114!! use netcdf, only: &
115!! & NF90_NOERR, NF90_ENOMEM, NF90_INQUIRE_VARIABLE
116 type(GD_NC_VARIABLE_ENTRY), intent(inout):: ent
117 integer, intent(out):: stat
118 integer:: ndims
119 if (ent%varid > 0) then
120 stat = nf90_inquire_variable(ent%fileid, ent%varid, ndims = ndims)
121 if (stat /= nf90_noerr) return
122 if ((ent%dimid > 0) .and. (ndims /= 1)) goto 100
123 if (ndims == 0) then
124 nullify(ent%dimids)
125 stat = nf90_noerr
126 return
127 endif
128 allocate(ent%dimids(ndims), stat=stat)
129 if (stat /= 0) then
130 stat = nf90_enomem
131 return
132 endif
133 stat = nf90_inquire_variable(ent%fileid, ent%varid, dimids = ent%dimids)
134 if (stat /= nf90_noerr) return
135 if ((ent%dimid > 0) .and. (ent%dimids(1) /= ent%dimid)) then
136 deallocate(ent%dimids)
137 goto 100
138 endif
139 else
140 allocate(ent%dimids(1), stat=stat)
141 if (stat /= 0) then
142 stat = nf90_enomem
143 return
144 endif
145 ent%dimids(1) = ent%dimid
146 endif
147 stat = nf90_noerr
148 return
149
150100 continue
151 ent%varid = 0
152 allocate(ent%dimids(1))
153 ent%dimids(1) = ent%dimid
154 end subroutine internal_build_dimids
155

◆ vtable_delete()

integer function, public gtdata_netcdf_internal::vtable_delete ( type(gd_nc_variable), intent(in)  var)

Definition at line 160 of file gtdata_netcdf_internal.f90.

161 type(GD_NC_VARIABLE), intent(in):: var
162 if (.not. allocated(gdnctab)) goto 999
163 if (var%id <= 0 .or. var%id > size(gdnctab)) goto 999
164 if (gdnctab(var%id)%fileid == 0) goto 999
165 result = gdnctab(var%id)%fileid
166 gdnctab(var%id)%fileid = 0
167 gdnctab(var%id)%varid = 0
168 gdnctab(var%id)%dimid = 0
169 gdnctab(var%id)%attrid = 0
170 if (associated(gdnctab(var%id)%dimids)) &
171 & deallocate(gdnctab(var%id)%dimids)
172 call dbgmessage('gtdata_netcdf_internal.delete: delete %d', i=(/var%id/))
173 return
174 !
175999 continue
176 result = nf90_enotvar

◆ vtable_lookup()

integer function, public gtdata_netcdf_internal::vtable_lookup ( type(gd_nc_variable), intent(in)  var,
type(gd_nc_variable_entry), intent(out)  entry 
)

Definition at line 179 of file gtdata_netcdf_internal.f90.

180 type(GD_NC_VARIABLE), intent(in):: var
181 type(GD_NC_VARIABLE_ENTRY), intent(out):: entry
182 if (.not. allocated(gdnctab)) goto 999
183 if (var%id <= 0 .or. var%id > size(gdnctab)) goto 999
184 if (gdnctab(var%id)%fileid == 0) goto 999
185 entry = gdnctab(var%id)
186 result = nf90_noerr
187 return
188 !
189999 continue
190 nullify(entry%dimids)
191 entry%fileid = -1
192 entry%varid = -1
193 entry%dimid = -1
194 entry%attrid = -1
195 result = nf90_enotvar

◆ vtable_set_attrid()

integer function, public gtdata_netcdf_internal::vtable_set_attrid ( type(gd_nc_variable), intent(in)  var,
integer, intent(in)  attrid 
)

Definition at line 198 of file gtdata_netcdf_internal.f90.

199 type(GD_NC_VARIABLE), intent(in):: var
200 integer, intent(in):: attrid
201 continue
202 if (.not. allocated(gdnctab)) goto 999
203 if (var%id <= 0 .or. var%id > size(gdnctab)) goto 999
204 if (gdnctab(var%id)%fileid == 0) goto 999
205 gdnctab(var%id)%attrid = attrid
206 result = nf90_noerr
207 return
208 !
209999 continue
210 result = nf90_enotvar