V3FIT
ezspline_save.f90
1 !/////
2 ! R8 !
3 !/////
4 !
5 ! MOD DMC Apr 2007 -- support Hybrid spline objects.
6 ! these objects can have "zonal step function" dimensions 1 less than
7 ! the grid dimension; also, the number of spline coefficients varies
8 ! depending on how many dimensions are using cubic interpolation. So,
9 ! local variables in0, in1, in2, in3 are defined to distinguish these
10 ! dimensions from the grid dimensions spline_o%n1, etc.
11 !
12 ! 1-D
13 !
14 
15 #if defined(NETCDF)
16  subroutine ezspline_save1_r8(spline_o, filename, ier, &
17  spl_name, fullsave)
18  use ezspline_obj
19  use ezcdf
20  implicit none
21  type(EZspline1_r8) :: spline_o
22  character*(*) :: filename
23  ! ier:
24  ! 17=could not save spline object in file filename
25  integer, intent(out) :: ier
26 
27  ! if SPL_NAME is set, APPEND to file instead of creating new; prepend
28  ! name to all NetCDF data items. This allows one file to contain
29  ! multiple items. (new dmc Mar 2006)
30  character*(*), intent(in), OPTIONAL :: SPL_NAME ! optional spline "name"
31 
32  ! if FULLSAVE is set .TRUE., save derived coefficients along with data
33  ! this saves recalculating the coefficients when file is read (Mar 2006)
34  logical, intent(in), OPTIONAL :: FULLSAVE
35 
36  integer ncid, ifail, in0, in1
37  integer dimlens(3)
38  character(2), parameter :: real8='R8'
39  character(3), parameter :: int='INT'
40 
41  character*32 zpre
42  logical :: fullsv,imodify
43 
44  ier = 0
45  if(spline_o%isReady /= 1) then
46  ier = 94
47  return
48  endif
49 
50  in0 = size(spline_o%fspl,1)
51  in1 = size(spline_o%fspl,2)
52 
53  if(present(spl_name)) then
54  call ezspline_spl_name_chk(spl_name,ier)
55  if(ier.ne.0) return
56  zpre=spl_name
57  else
58  zpre=' '
59  endif
60 
61  if(present(fullsave)) then
62  fullsv=fullsave
63  else
64  fullsv=.false.
65  endif
66 
67  if(zpre.eq.' ') then
68  call cdfopn(ncid, filename, 'w') ! no error flag??
69  imodify=.false.
70  else
71  call cdfopn(ncid, filename, 'v') ! Append if exists, else create new
72  imodify=.true.
73  endif
74  if(ncid==0) then
75  ier = 39
76  return
77  endif
78 
79  dimlens = (/1, 1, 1/) ! scalar
80  call ezspline_defvar(ncid, trim(zpre)//'n1', dimlens, int, imodify, ier)
81  call ezspline_defvar(ncid, trim(zpre)//'klookup1', dimlens, int, imodify, ier)
82  call ezspline_defvar(ncid, trim(zpre)//'isHermite', dimlens, int, imodify, ier)
83  call ezspline_defvar(ncid, trim(zpre)//'isLinear', dimlens, int, imodify, ier)
84  call ezspline_defvar(ncid, trim(zpre)//'isReady', dimlens, int, imodify, ier)
85  dimlens = (/2, 1, 1/)
86  call ezspline_defvar(ncid, trim(zpre)//'ibctype1', dimlens, int, imodify, ier)
87  dimlens = (/spline_o%n1, 1, 1/)
88  call ezspline_defvar(ncid, trim(zpre)//'x1', dimlens, real8, imodify, ier)
89  dimlens = (/1, 1, 1/)
90  call ezspline_defvar(ncid, trim(zpre)//'bcval1min', dimlens, real8, imodify, ier)
91  call ezspline_defvar(ncid, trim(zpre)//'bcval1max', dimlens, real8, imodify, ier)
92  !
93  ! save only the function values at the nodes
94  !
95  dimlens = (/spline_o%n1, 1, 1/)
96  if(.not.fullsv) then
97  call ezspline_defvar(ncid, trim(zpre)//'f', dimlens, real8, imodify, ier)
98  else
99  dimlens = (/1, 1, 1/) ! scalar
100  call ezspline_defvar(ncid, trim(zpre)//'x1min', dimlens, real8, imodify, ier)
101  call ezspline_defvar(ncid, trim(zpre)//'x1max', dimlens, real8, imodify, ier)
102  call ezspline_defvar(ncid, trim(zpre)//'ilin1', dimlens, int, imodify, ier)
103  dimlens = (/spline_o%n1, 4, 1/)
104  call ezspline_defvar(ncid, trim(zpre)//'x1pkg', dimlens, real8, imodify, ier)
105  dimlens = (/2, spline_o%n1, 1/)
106  call ezspline_defvar(ncid, trim(zpre)//'fspl', dimlens, real8, imodify, ier)
107  endif
108  if(ier.ne.0) return
109 
110  call cdfputvar(ncid, trim(zpre)//'n1', spline_o%n1, ifail)
111  call cdfputvar(ncid, trim(zpre)//'klookup1', spline_o%klookup1, ifail)
112  call cdfputvar(ncid, trim(zpre)//'isHermite', spline_o%isHermite, ifail)
113  call cdfputvar(ncid, trim(zpre)//'isLinear', spline_o%isLinear, ifail)
114  call cdfputvar(ncid, trim(zpre)//'isReady', spline_o%isReady, ifail)
115  call cdfputvar(ncid, trim(zpre)//'ibctype1', spline_o%ibctype1, ifail)
116  call cdfputvar(ncid, trim(zpre)//'x1', spline_o%x1, ifail)
117  call cdfputvar(ncid, trim(zpre)//'bcval1min', spline_o%bcval1min, ifail)
118  call cdfputvar(ncid, trim(zpre)//'bcval1max', spline_o%bcval1max, ifail)
119  if(ifail/=0) then
120  ier=40
121  return
122  endif
123 
124  if(spline_o%isReady == 1) then
125  if(.not.fullsv) then
126  call cdfputvar(ncid, trim(zpre)//'f', spline_o%fspl(1,:), ifail)
127  else
128  call cdfputvar(ncid, trim(zpre)//'x1min', spline_o%x1min, ifail)
129  call cdfputvar(ncid, trim(zpre)//'x1max', spline_o%x1max, ifail)
130  call cdfputvar(ncid, trim(zpre)//'ilin1', spline_o%ilin1, ifail)
131  call cdfputvar(ncid, trim(zpre)//'x1pkg', spline_o%x1pkg, ifail)
132  call cdfputvar(ncid, trim(zpre)//'fspl', spline_o%fspl, ifail)
133  endif
134  else
135  ier = 93
136  endif
137 
138  call cdfcls(ncid)
139 
140  if(ifail /=0) ier = 17
141 
142  return
143  end subroutine ezspline_save1_r8
144 
145 !!
146 !! 2-D
147 !!
148 
149 
150  subroutine ezspline_save2_r8(spline_o, filename, ier, &
151  spl_name, fullsave)
152  use ezspline_obj
153  use ezcdf
154  implicit none
155  type(EZspline2_r8) :: spline_o
156  character*(*) :: filename
157  ! ier:
158  ! 17=could not save spline object in file filename
159  integer, intent(out) :: ier
160 
161  ! if SPL_NAME is set, APPEND to file instead of creating new; prepend
162  ! name to all NetCDF data items. This allows one file to contain
163  ! multiple items. (new dmc Mar 2006)
164  character*(*), intent(in), OPTIONAL :: SPL_NAME ! optional spline "name"
165 
166  ! if FULLSAVE is set .TRUE., save derived coefficients along with data
167  ! this saves recalculating the coefficients when file is read (Mar 2006)
168  logical, intent(in), OPTIONAL :: FULLSAVE
169 
170  integer ncid, ifail, in0, in1, in2
171  integer dimlens(3)
172  character(2), parameter :: real8='R8'
173  character(3), parameter :: int='INT'
174 
175  character*32 zpre
176  logical :: fullsv,imodify
177 
178  ier = 0
179  if(spline_o%isReady /= 1) then
180  ier = 94
181  return
182  endif
183 
184  in0 = size(spline_o%fspl,1)
185  in1 = size(spline_o%fspl,2)
186  in2 = size(spline_o%fspl,3)
187 
188  if(present(spl_name)) then
189  call ezspline_spl_name_chk(spl_name,ier)
190  if(ier.ne.0) return
191  zpre=spl_name
192  else
193  zpre=' '
194  endif
195 
196  if(present(fullsave)) then
197  fullsv=fullsave
198  else
199  fullsv=.false.
200  endif
201 
202  if(zpre.eq.' ') then
203  call cdfopn(ncid, filename, 'w') ! no error flag??
204  imodify=.false.
205  else
206  call cdfopn(ncid, filename, 'v') ! Append if exists, else create new
207  imodify=.true.
208  endif
209  if(ncid==0) then
210  ier = 39
211  return
212  endif
213 
214  dimlens = (/1, 1, 1/) ! scalar
215  call ezspline_defvar(ncid, trim(zpre)//'n1', dimlens, int, imodify, ier)
216  call ezspline_defvar(ncid, trim(zpre)//'n2', dimlens, int, imodify, ier)
217  call ezspline_defvar(ncid, trim(zpre)//'klookup1', dimlens, int, imodify, ier)
218  call ezspline_defvar(ncid, trim(zpre)//'klookup2', dimlens, int, imodify, ier)
219  call ezspline_defvar(ncid, trim(zpre)//'isHermite', dimlens, int, imodify, ier)
220  call ezspline_defvar(ncid, trim(zpre)//'isLinear', dimlens, int, imodify, ier)
221  call ezspline_defvar(ncid, trim(zpre)//'isHybrid', dimlens, int, imodify, ier)
222  call ezspline_defvar(ncid, trim(zpre)//'isReady', dimlens, int, imodify, ier)
223  dimlens = (/2, 1, 1/)
224  call ezspline_defvar(ncid, trim(zpre)//'hspline', dimlens, int, imodify, ier)
225  dimlens = (/2, 1, 1/)
226  call ezspline_defvar(ncid, trim(zpre)//'ibctype1', dimlens, int, imodify, ier)
227  call ezspline_defvar(ncid, trim(zpre)//'ibctype2', dimlens, int, imodify, ier)
228  dimlens = (/spline_o%n1, 1, 1/)
229  call ezspline_defvar(ncid, trim(zpre)//'x1', dimlens, real8, imodify, ier)
230  dimlens = (/spline_o%n2, 1, 1/)
231  call ezspline_defvar(ncid, trim(zpre)//'x2', dimlens, real8, imodify, ier)
232  dimlens = (/in2, 1, 1/)
233  call ezspline_defvar(ncid, trim(zpre)//'bcval1min', dimlens, real8, imodify, ier)
234  call ezspline_defvar(ncid, trim(zpre)//'bcval1max', dimlens, real8, imodify, ier)
235  dimlens = (/in1, 1, 1/)
236  call ezspline_defvar(ncid, trim(zpre)//'bcval2min', dimlens, real8, imodify, ier)
237  call ezspline_defvar(ncid, trim(zpre)//'bcval2max', dimlens, real8, imodify, ier)
238  !
239  ! save only the function values at the nodes
240  !
241  dimlens = (/in1, in2, 1/)
242  if(.not.fullsv) then
243  call ezspline_defvar(ncid, trim(zpre)//'f', dimlens, real8, imodify, ier)
244  else
245  dimlens = (/1, 1, 1/) ! scalar
246  call ezspline_defvar(ncid, trim(zpre)//'x1min', dimlens, real8, imodify, ier)
247  call ezspline_defvar(ncid, trim(zpre)//'x1max', dimlens, real8, imodify, ier)
248  call ezspline_defvar(ncid, trim(zpre)//'ilin1', dimlens, int, imodify, ier)
249  call ezspline_defvar(ncid, trim(zpre)//'x2min', dimlens, real8, imodify, ier)
250  call ezspline_defvar(ncid, trim(zpre)//'x2max', dimlens, real8, imodify, ier)
251  call ezspline_defvar(ncid, trim(zpre)//'ilin2', dimlens, int, imodify, ier)
252  dimlens = (/spline_o%n1, 4, 1/)
253  call ezspline_defvar(ncid, trim(zpre)//'x1pkg', dimlens, real8, imodify, ier)
254  dimlens = (/spline_o%n2, 4, 1/)
255  call ezspline_defvar(ncid, trim(zpre)//'x2pkg', dimlens, real8, imodify, ier)
256  dimlens = (/in0, in1, in2/)
257  call ezspline_defvar(ncid, trim(zpre)//'fspl', dimlens, real8, imodify, ier)
258  endif
259  if(ier.ne.0) return
260 
261  call cdfputvar(ncid, trim(zpre)//'n1', spline_o%n1, ifail)
262  call cdfputvar(ncid, trim(zpre)//'n2', spline_o%n2, ifail)
263  call cdfputvar(ncid, trim(zpre)//'klookup1', spline_o%klookup1, ifail)
264  call cdfputvar(ncid, trim(zpre)//'klookup2', spline_o%klookup2, ifail)
265  call cdfputvar(ncid, trim(zpre)//'isHermite', spline_o%isHermite, ifail)
266  call cdfputvar(ncid, trim(zpre)//'isLinear', spline_o%isLinear, ifail)
267  call cdfputvar(ncid, trim(zpre)//'isHybrid', spline_o%isHybrid, ifail)
268  call cdfputvar(ncid, trim(zpre)//'isReady', spline_o%isReady, ifail)
269  call cdfputvar(ncid, trim(zpre)//'hspline', spline_o%hspline, ifail)
270  call cdfputvar(ncid, trim(zpre)//'ibctype1', spline_o%ibctype1, ifail)
271  call cdfputvar(ncid, trim(zpre)//'ibctype2', spline_o%ibctype2, ifail)
272  call cdfputvar(ncid, trim(zpre)//'x1', spline_o%x1, ifail)
273  call cdfputvar(ncid, trim(zpre)//'x2', spline_o%x2, ifail)
274  call cdfputvar(ncid, trim(zpre)//'bcval1min', spline_o%bcval1min, ifail)
275  call cdfputvar(ncid, trim(zpre)//'bcval1max', spline_o%bcval1max, ifail)
276  call cdfputvar(ncid, trim(zpre)//'bcval2min', spline_o%bcval2min, ifail)
277  call cdfputvar(ncid, trim(zpre)//'bcval2max', spline_o%bcval2max, ifail)
278  if(ifail/=0) then
279  ier=40
280  return
281  endif
282 
283  if(spline_o%isReady == 1) then
284  if(.not.fullsv) then
285  call cdfputvar(ncid, trim(zpre)//'f', spline_o%fspl(1,:,:), ifail)
286  else
287  call cdfputvar(ncid, trim(zpre)//'x1min', spline_o%x1min, ifail)
288  call cdfputvar(ncid, trim(zpre)//'x1max', spline_o%x1max, ifail)
289  call cdfputvar(ncid, trim(zpre)//'ilin1', spline_o%ilin1, ifail)
290  call cdfputvar(ncid, trim(zpre)//'x2min', spline_o%x2min, ifail)
291  call cdfputvar(ncid, trim(zpre)//'x2max', spline_o%x2max, ifail)
292  call cdfputvar(ncid, trim(zpre)//'ilin2', spline_o%ilin2, ifail)
293  call cdfputvar(ncid, trim(zpre)//'x1pkg', spline_o%x1pkg, ifail)
294  call cdfputvar(ncid, trim(zpre)//'x2pkg', spline_o%x2pkg, ifail)
295  call cdfputvar(ncid, trim(zpre)//'fspl', spline_o%fspl, ifail)
296  endif
297  else
298  ier = 93
299  endif
300 
301  call cdfcls(ncid)
302 
303  if(ifail /=0) ier = 17
304 
305  return
306  end subroutine ezspline_save2_r8
307 
308 
309 
310 
311 
312 !!!
313 !!! 3-D
314 !!!
315 
316 
317 
318  subroutine ezspline_save3_r8(spline_o, filename, ier, &
319  spl_name, fullsave)
320  use ezspline_obj
321  use ezcdf
322  implicit none
323  type(EZspline3_r8) :: spline_o
324  character*(*) :: filename
325  ! ier:
326  ! 17=could not save spline object in file filename
327  integer, intent(out) :: ier
328 
329  ! if SPL_NAME is set, APPEND to file instead of creating new; prepend
330  ! name to all NetCDF data items. This allows one file to contain
331  ! multiple items. (new dmc Mar 2006)
332  character*(*), intent(in), OPTIONAL :: SPL_NAME ! optional spline "name"
333 
334  ! if FULLSAVE is set .TRUE., save derived coefficients along with data
335  ! this saves recalculating the coefficients when file is read (Mar 2006)
336  logical, intent(in), OPTIONAL :: FULLSAVE
337 
338  integer ncid, ifail, in0, in1, in2, in3
339  integer dimlens(3)
340  character(2), parameter :: real8='R8'
341  character(3), parameter :: int='INT'
342 
343  character*32 zpre
344  logical :: fullsv,imodify
345 
346  ier = 0
347  if(spline_o%isReady /= 1) then
348  ier = 94
349  return
350  endif
351 
352  in0 = size(spline_o%fspl,1)
353  in1 = size(spline_o%fspl,2)
354  in2 = size(spline_o%fspl,3)
355  in3 = size(spline_o%fspl,4)
356 
357  if(present(spl_name)) then
358  call ezspline_spl_name_chk(spl_name,ier)
359  if(ier.ne.0) return
360  zpre=spl_name
361  else
362  zpre=' '
363  endif
364 
365  if(present(fullsave)) then
366  fullsv=fullsave
367  else
368  fullsv=.false.
369  endif
370 
371  if(zpre.eq.' ') then
372  call cdfopn(ncid, filename, 'w') ! no error flag??
373  imodify=.false.
374  else
375  call cdfopn(ncid, filename, 'v') ! Append if exists, else create new
376  imodify=.true.
377  endif
378  if(ncid==0) then
379  ier = 39
380  return
381  endif
382 
383  dimlens = (/1, 1, 1/) ! scalar
384  call ezspline_defvar(ncid, trim(zpre)//'n1', dimlens, int, imodify, ier)
385  call ezspline_defvar(ncid, trim(zpre)//'n2', dimlens, int, imodify, ier)
386  call ezspline_defvar(ncid, trim(zpre)//'n3', dimlens, int, imodify, ier)
387  call ezspline_defvar(ncid, trim(zpre)//'klookup1', dimlens, int, imodify, ier)
388  call ezspline_defvar(ncid, trim(zpre)//'klookup2', dimlens, int, imodify, ier)
389  call ezspline_defvar(ncid, trim(zpre)//'klookup3', dimlens, int, imodify, ier)
390  call ezspline_defvar(ncid, trim(zpre)//'isHermite', dimlens, int, imodify, ier)
391  call ezspline_defvar(ncid, trim(zpre)//'isLinear', dimlens, int, imodify, ier)
392  call ezspline_defvar(ncid, trim(zpre)//'isHybrid', dimlens, int, imodify, ier)
393  call ezspline_defvar(ncid, trim(zpre)//'isReady', dimlens, int, imodify, ier)
394  dimlens = (/3, 1, 1/)
395  call ezspline_defvar(ncid, trim(zpre)//'hspline', dimlens, int, imodify, ier)
396  dimlens = (/2, 1, 1/)
397  call ezspline_defvar(ncid, trim(zpre)//'ibctype1', dimlens, int, imodify, ier)
398  call ezspline_defvar(ncid, trim(zpre)//'ibctype2', dimlens, int, imodify, ier)
399  call ezspline_defvar(ncid, trim(zpre)//'ibctype3', dimlens, int, imodify, ier)
400  dimlens = (/spline_o%n1, 1, 1/)
401  call ezspline_defvar(ncid, trim(zpre)//'x1', dimlens, real8, imodify, ier)
402  dimlens = (/spline_o%n2, 1, 1/)
403  call ezspline_defvar(ncid, trim(zpre)//'x2', dimlens, real8, imodify, ier)
404  dimlens = (/spline_o%n3, 1, 1/)
405  call ezspline_defvar(ncid, trim(zpre)//'x3', dimlens, real8, imodify, ier)
406  dimlens = (/in2, in3, 1/)
407  call ezspline_defvar(ncid, trim(zpre)//'bcval1min', dimlens, real8, imodify, ier)
408  call ezspline_defvar(ncid, trim(zpre)//'bcval1max', dimlens, real8, imodify, ier)
409  dimlens = (/in1, in3, 1/)
410  call ezspline_defvar(ncid, trim(zpre)//'bcval2min', dimlens, real8, imodify, ier)
411  call ezspline_defvar(ncid, trim(zpre)//'bcval2max', dimlens, real8, imodify, ier)
412  dimlens = (/in1, in2, 1/)
413  call ezspline_defvar(ncid, trim(zpre)//'bcval3min', dimlens, real8, imodify, ier)
414  call ezspline_defvar(ncid, trim(zpre)//'bcval3max', dimlens, real8, imodify, ier)
415  !
416  ! save only the function values at the nodes
417  !
418  dimlens = (/in1, in2, in3/)
419  if(.not.fullsv) then
420  call ezspline_defvar(ncid, trim(zpre)//'f', dimlens, real8, imodify, ier)
421  else
422  dimlens = (/1, 1, 1/) ! scalar
423  call ezspline_defvar(ncid, trim(zpre)//'x1min', dimlens, real8, imodify, ier)
424  call ezspline_defvar(ncid, trim(zpre)//'x1max', dimlens, real8, imodify, ier)
425  call ezspline_defvar(ncid, trim(zpre)//'ilin1', dimlens, int, imodify, ier)
426  call ezspline_defvar(ncid, trim(zpre)//'x2min', dimlens, real8, imodify, ier)
427  call ezspline_defvar(ncid, trim(zpre)//'x2max', dimlens, real8, imodify, ier)
428  call ezspline_defvar(ncid, trim(zpre)//'ilin2', dimlens, int, imodify, ier)
429  call ezspline_defvar(ncid, trim(zpre)//'x3min', dimlens, real8, imodify, ier)
430  call ezspline_defvar(ncid, trim(zpre)//'x3max', dimlens, real8, imodify, ier)
431  call ezspline_defvar(ncid, trim(zpre)//'ilin3', dimlens, int, imodify, ier)
432  dimlens = (/spline_o%n1, 4, 1/)
433  call ezspline_defvar(ncid, trim(zpre)//'x1pkg', dimlens, real8, imodify, ier)
434  dimlens = (/spline_o%n2, 4, 1/)
435  call ezspline_defvar(ncid, trim(zpre)//'x2pkg', dimlens, real8, imodify, ier)
436  dimlens = (/spline_o%n3, 4, 1/)
437  call ezspline_defvar(ncid, trim(zpre)//'x3pkg', dimlens, real8, imodify, ier)
438  dimlens = (/in0*in1, in2, in3/)
439  call ezspline_defvar(ncid, trim(zpre)//'fspl', dimlens, real8, imodify, ier)
440  endif
441  if(ier.ne.0) return
442 
443  call cdfputvar(ncid, trim(zpre)//'n1', spline_o%n1, ifail)
444  call cdfputvar(ncid, trim(zpre)//'n2', spline_o%n2, ifail)
445  call cdfputvar(ncid, trim(zpre)//'n3', spline_o%n3, ifail)
446  call cdfputvar(ncid, trim(zpre)//'klookup1', spline_o%klookup1, ifail)
447  call cdfputvar(ncid, trim(zpre)//'klookup2', spline_o%klookup2, ifail)
448  call cdfputvar(ncid, trim(zpre)//'klookup3', spline_o%klookup3, ifail)
449  call cdfputvar(ncid, trim(zpre)//'isHermite', spline_o%isHermite, ifail)
450  call cdfputvar(ncid, trim(zpre)//'isLinear', spline_o%isLinear, ifail)
451  call cdfputvar(ncid, trim(zpre)//'isHybrid', spline_o%isHybrid, ifail)
452  call cdfputvar(ncid, trim(zpre)//'isReady', spline_o%isReady, ifail)
453  call cdfputvar(ncid, trim(zpre)//'hspline', spline_o%hspline, ifail)
454  call cdfputvar(ncid, trim(zpre)//'ibctype1', spline_o%ibctype1, ifail)
455  call cdfputvar(ncid, trim(zpre)//'ibctype2', spline_o%ibctype2, ifail)
456  call cdfputvar(ncid, trim(zpre)//'ibctype3', spline_o%ibctype3, ifail)
457  call cdfputvar(ncid, trim(zpre)//'x1', spline_o%x1, ifail)
458  call cdfputvar(ncid, trim(zpre)//'x2', spline_o%x2, ifail)
459  call cdfputvar(ncid, trim(zpre)//'x3', spline_o%x3, ifail)
460  call cdfputvar(ncid, trim(zpre)//'bcval1min', spline_o%bcval1min, ifail)
461  call cdfputvar(ncid, trim(zpre)//'bcval1max', spline_o%bcval1max, ifail)
462  call cdfputvar(ncid, trim(zpre)//'bcval2min', spline_o%bcval2min, ifail)
463  call cdfputvar(ncid, trim(zpre)//'bcval2max', spline_o%bcval2max, ifail)
464  call cdfputvar(ncid, trim(zpre)//'bcval3min', spline_o%bcval3min, ifail)
465  call cdfputvar(ncid, trim(zpre)//'bcval3max', spline_o%bcval3max, ifail)
466  if(ifail/=0) then
467  ier=40
468  return
469  endif
470 
471  if(spline_o%isReady == 1) then
472  if(.not.fullsv) then
473  call cdfputvar(ncid, trim(zpre)//'f', spline_o%fspl(1,:,:,:), ifail)
474  else
475  call cdfputvar(ncid, trim(zpre)//'x1min', spline_o%x1min, ifail)
476  call cdfputvar(ncid, trim(zpre)//'x1max', spline_o%x1max, ifail)
477  call cdfputvar(ncid, trim(zpre)//'ilin1', spline_o%ilin1, ifail)
478  call cdfputvar(ncid, trim(zpre)//'x2min', spline_o%x2min, ifail)
479  call cdfputvar(ncid, trim(zpre)//'x2max', spline_o%x2max, ifail)
480  call cdfputvar(ncid, trim(zpre)//'ilin2', spline_o%ilin2, ifail)
481  call cdfputvar(ncid, trim(zpre)//'x3min', spline_o%x3min, ifail)
482  call cdfputvar(ncid, trim(zpre)//'x3max', spline_o%x3max, ifail)
483  call cdfputvar(ncid, trim(zpre)//'ilin3', spline_o%ilin3, ifail)
484  call cdfputvar(ncid, trim(zpre)//'x1pkg', spline_o%x1pkg, ifail)
485  call cdfputvar(ncid, trim(zpre)//'x2pkg', spline_o%x2pkg, ifail)
486  call cdfputvar(ncid, trim(zpre)//'x3pkg', spline_o%x3pkg, ifail)
487  call ezspline_cdfput3(ncid, trim(zpre)//'fspl', spline_o%fspl, &
488  in0*in1, in2, in3, ifail)
489  endif
490  else
491  ier = 93
492  endif
493 
494  call cdfcls(ncid)
495 
496  if(ifail /=0) ier = 17
497 
498  return
499  end subroutine ezspline_save3_r8
500 !/////
501 ! R4 !
502 !/////
503 !
504 ! 1-D
505 !
506 
507 
508  subroutine ezspline_save1_r4(spline_o, filename, ier, &
509  spl_name, fullsave)
510  use ezspline_obj
511  use ezcdf
512  implicit none
513  type(EZspline1_r4) :: spline_o
514  character*(*) :: filename
515  ! ier:
516  ! 17=could not save spline object in file filename
517  integer, intent(out) :: ier
518 
519  ! if SPL_NAME is set, APPEND to file instead of creating new; prepend
520  ! name to all NetCDF data items. This allows one file to contain
521  ! multiple items. (new dmc Mar 2006)
522  character*(*), intent(in), OPTIONAL :: SPL_NAME ! optional spline "name"
523 
524  ! if FULLSAVE is set .TRUE., save derived coefficients along with data
525  ! this saves recalculating the coefficients when file is read (Mar 2006)
526  logical, intent(in), OPTIONAL :: FULLSAVE
527 
528  integer ncid, ifail, in0, in1
529  integer dimlens(3)
530  character(2), parameter :: real4='R4'
531  character(3), parameter :: int='INT'
532 
533  character*32 zpre
534  logical :: fullsv,imodify
535 
536  ier = 0
537  if(spline_o%isReady /= 1) then
538  ier = 94
539  return
540  endif
541 
542  in0 = size(spline_o%fspl,1)
543  in1 = size(spline_o%fspl,2)
544 
545  if(present(spl_name)) then
546  call ezspline_spl_name_chk(spl_name,ier)
547  if(ier.ne.0) return
548  zpre=spl_name
549  else
550  zpre=' '
551  endif
552 
553  if(present(fullsave)) then
554  fullsv=fullsave
555  else
556  fullsv=.false.
557  endif
558 
559  if(zpre.eq.' ') then
560  call cdfopn(ncid, filename, 'w') ! no error flag??
561  imodify=.false.
562  else
563  call cdfopn(ncid, filename, 'v') ! Append if exists, else create new
564  imodify=.true.
565  endif
566  if(ncid==0) then
567  ier = 39
568  return
569  endif
570 
571  dimlens = (/1, 1, 1/) ! scalar
572  call ezspline_defvar(ncid, trim(zpre)//'n1', dimlens, int, imodify, ier)
573  call ezspline_defvar(ncid, trim(zpre)//'klookup1', dimlens, int, imodify, ier)
574  call ezspline_defvar(ncid, trim(zpre)//'isHermite', dimlens, int, imodify, ier)
575  call ezspline_defvar(ncid, trim(zpre)//'isLinear', dimlens, int, imodify, ier)
576  call ezspline_defvar(ncid, trim(zpre)//'isReady', dimlens, int, imodify, ier)
577  dimlens = (/2, 1, 1/)
578  call ezspline_defvar(ncid, trim(zpre)//'ibctype1', dimlens, int, imodify, ier)
579  dimlens = (/spline_o%n1, 1, 1/)
580  call ezspline_defvar(ncid, trim(zpre)//'x1', dimlens, real4, imodify, ier)
581  dimlens = (/1, 1, 1/)
582  call ezspline_defvar(ncid, trim(zpre)//'bcval1min', dimlens, real4, imodify, ier)
583  call ezspline_defvar(ncid, trim(zpre)//'bcval1max', dimlens, real4, imodify, ier)
584  !
585  ! save only the function values at the nodes
586  !
587  dimlens = (/spline_o%n1, 1, 1/)
588  if(.not.fullsv) then
589  call ezspline_defvar(ncid, trim(zpre)//'f', dimlens, real4, imodify, ier)
590  else
591  dimlens = (/1, 1, 1/) ! scalar
592  call ezspline_defvar(ncid, trim(zpre)//'x1min', dimlens, real4, imodify, ier)
593  call ezspline_defvar(ncid, trim(zpre)//'x1max', dimlens, real4, imodify, ier)
594  call ezspline_defvar(ncid, trim(zpre)//'ilin1', dimlens, int, imodify, ier)
595  dimlens = (/spline_o%n1, 4, 1/)
596  call ezspline_defvar(ncid, trim(zpre)//'x1pkg', dimlens, real4, imodify, ier)
597  dimlens = (/2, spline_o%n1, 1/)
598  call ezspline_defvar(ncid, trim(zpre)//'fspl', dimlens, real4, imodify, ier)
599  endif
600  if(ier.ne.0) return
601 
602  call cdfputvar(ncid, trim(zpre)//'n1', spline_o%n1, ifail)
603  call cdfputvar(ncid, trim(zpre)//'klookup1', spline_o%klookup1, ifail)
604  call cdfputvar(ncid, trim(zpre)//'isHermite', spline_o%isHermite, ifail)
605  call cdfputvar(ncid, trim(zpre)//'isLinear', spline_o%isLinear, ifail)
606  call cdfputvar(ncid, trim(zpre)//'isReady', spline_o%isReady, ifail)
607  call cdfputvar(ncid, trim(zpre)//'ibctype1', spline_o%ibctype1, ifail)
608  call cdfputvar(ncid, trim(zpre)//'x1', spline_o%x1, ifail)
609  call cdfputvar(ncid, trim(zpre)//'bcval1min', spline_o%bcval1min, ifail)
610  call cdfputvar(ncid, trim(zpre)//'bcval1max', spline_o%bcval1max, ifail)
611  if(ifail/=0) then
612  ier=40
613  return
614  endif
615 
616  if(spline_o%isReady == 1) then
617  if(.not.fullsv) then
618  call cdfputvar(ncid, trim(zpre)//'f', spline_o%fspl(1,:), ifail)
619  else
620  call cdfputvar(ncid, trim(zpre)//'x1min', spline_o%x1min, ifail)
621  call cdfputvar(ncid, trim(zpre)//'x1max', spline_o%x1max, ifail)
622  call cdfputvar(ncid, trim(zpre)//'ilin1', spline_o%ilin1, ifail)
623  call cdfputvar(ncid, trim(zpre)//'x1pkg', spline_o%x1pkg, ifail)
624  call cdfputvar(ncid, trim(zpre)//'fspl', spline_o%fspl, ifail)
625  endif
626  else
627  ier = 93
628  endif
629 
630  call cdfcls(ncid)
631 
632  if(ifail /=0) ier = 17
633 
634  return
635  end subroutine ezspline_save1_r4
636 
637 !!
638 !! 2-D
639 !!
640 
641 
642  subroutine ezspline_save2_r4(spline_o, filename, ier, &
643  spl_name, fullsave)
644  use ezspline_obj
645  use ezcdf
646  implicit none
647  type(EZspline2_r4) :: spline_o
648  character*(*) :: filename
649  ! ier:
650  ! 17=could not save spline object in file filename
651  integer, intent(out) :: ier
652 
653  ! if SPL_NAME is set, APPEND to file instead of creating new; prepend
654  ! name to all NetCDF data items. This allows one file to contain
655  ! multiple items. (new dmc Mar 2006)
656  character*(*), intent(in), OPTIONAL :: SPL_NAME ! optional spline "name"
657 
658  ! if FULLSAVE is set .TRUE., save derived coefficients along with data
659  ! this saves recalculating the coefficients when file is read (Mar 2006)
660  logical, intent(in), OPTIONAL :: FULLSAVE
661 
662  integer ncid, ifail, in0, in1, in2
663  integer dimlens(3)
664  character(2), parameter :: real4='R4'
665  character(3), parameter :: int='INT'
666 
667  character*32 zpre
668  logical :: fullsv,imodify
669 
670  ier = 0
671  if(spline_o%isReady /= 1) then
672  ier = 94
673  return
674  endif
675 
676  in0 = size(spline_o%fspl,1)
677  in1 = size(spline_o%fspl,2)
678  in2 = size(spline_o%fspl,3)
679 
680  if(present(spl_name)) then
681  call ezspline_spl_name_chk(spl_name,ier)
682  if(ier.ne.0) return
683  zpre=spl_name
684  else
685  zpre=' '
686  endif
687 
688  if(present(fullsave)) then
689  fullsv=fullsave
690  else
691  fullsv=.false.
692  endif
693 
694  if(zpre.eq.' ') then
695  call cdfopn(ncid, filename, 'w') ! no error flag??
696  imodify=.false.
697  else
698  call cdfopn(ncid, filename, 'v') ! Append if exists, else create new
699  imodify=.true.
700  endif
701  if(ncid==0) then
702  ier = 39
703  return
704  endif
705 
706  dimlens = (/1, 1, 1/) ! scalar
707  call ezspline_defvar(ncid, trim(zpre)//'n1', dimlens, int, imodify, ier)
708  call ezspline_defvar(ncid, trim(zpre)//'n2', dimlens, int, imodify, ier)
709  call ezspline_defvar(ncid, trim(zpre)//'klookup1', dimlens, int, imodify, ier)
710  call ezspline_defvar(ncid, trim(zpre)//'klookup2', dimlens, int, imodify, ier)
711  call ezspline_defvar(ncid, trim(zpre)//'isHermite', dimlens, int, imodify, ier)
712  call ezspline_defvar(ncid, trim(zpre)//'isLinear', dimlens, int, imodify, ier)
713  call ezspline_defvar(ncid, trim(zpre)//'isHybrid', dimlens, int, imodify, ier)
714  call ezspline_defvar(ncid, trim(zpre)//'isReady', dimlens, int, imodify, ier)
715  dimlens = (/2, 1, 1/)
716  call ezspline_defvar(ncid, trim(zpre)//'hspline', dimlens, int, imodify, ier)
717  dimlens = (/2, 1, 1/)
718  call ezspline_defvar(ncid, trim(zpre)//'ibctype1', dimlens, int, imodify, ier)
719  call ezspline_defvar(ncid, trim(zpre)//'ibctype2', dimlens, int, imodify, ier)
720  dimlens = (/spline_o%n1, 1, 1/)
721  call ezspline_defvar(ncid, trim(zpre)//'x1', dimlens, real4, imodify, ier)
722  dimlens = (/spline_o%n2, 1, 1/)
723  call ezspline_defvar(ncid, trim(zpre)//'x2', dimlens, real4, imodify, ier)
724  dimlens = (/in2, 1, 1/)
725  call ezspline_defvar(ncid, trim(zpre)//'bcval1min', dimlens, real4, imodify, ier)
726  call ezspline_defvar(ncid, trim(zpre)//'bcval1max', dimlens, real4, imodify, ier)
727  dimlens = (/in1, 1, 1/)
728  call ezspline_defvar(ncid, trim(zpre)//'bcval2min', dimlens, real4, imodify, ier)
729  call ezspline_defvar(ncid, trim(zpre)//'bcval2max', dimlens, real4, imodify, ier)
730  !
731  ! save only the function values at the nodes
732  !
733  dimlens = (/in1, in2, 1/)
734  if(.not.fullsv) then
735  call ezspline_defvar(ncid, trim(zpre)//'f', dimlens, real4, imodify, ier)
736  else
737  dimlens = (/1, 1, 1/) ! scalar
738  call ezspline_defvar(ncid, trim(zpre)//'x1min', dimlens, real4, imodify, ier)
739  call ezspline_defvar(ncid, trim(zpre)//'x1max', dimlens, real4, imodify, ier)
740  call ezspline_defvar(ncid, trim(zpre)//'ilin1', dimlens, int, imodify, ier)
741  call ezspline_defvar(ncid, trim(zpre)//'x2min', dimlens, real4, imodify, ier)
742  call ezspline_defvar(ncid, trim(zpre)//'x2max', dimlens, real4, imodify, ier)
743  call ezspline_defvar(ncid, trim(zpre)//'ilin2', dimlens, int, imodify, ier)
744  dimlens = (/spline_o%n1, 4, 1/)
745  call ezspline_defvar(ncid, trim(zpre)//'x1pkg', dimlens, real4, imodify, ier)
746  dimlens = (/spline_o%n2, 4, 1/)
747  call ezspline_defvar(ncid, trim(zpre)//'x2pkg', dimlens, real4, imodify, ier)
748  dimlens = (/in0, in1, in2/)
749  call ezspline_defvar(ncid, trim(zpre)//'fspl', dimlens, real4, imodify, ier)
750  endif
751  if(ier.ne.0) return
752 
753  call cdfputvar(ncid, trim(zpre)//'n1', spline_o%n1, ifail)
754  call cdfputvar(ncid, trim(zpre)//'n2', spline_o%n2, ifail)
755  call cdfputvar(ncid, trim(zpre)//'klookup1', spline_o%klookup1, ifail)
756  call cdfputvar(ncid, trim(zpre)//'klookup2', spline_o%klookup2, ifail)
757  call cdfputvar(ncid, trim(zpre)//'isHermite', spline_o%isHermite, ifail)
758  call cdfputvar(ncid, trim(zpre)//'isLinear', spline_o%isLinear, ifail)
759  call cdfputvar(ncid, trim(zpre)//'isHybrid', spline_o%isHybrid, ifail)
760  call cdfputvar(ncid, trim(zpre)//'isReady', spline_o%isReady, ifail)
761  call cdfputvar(ncid, trim(zpre)//'hspline', spline_o%hspline, ifail)
762  call cdfputvar(ncid, trim(zpre)//'ibctype1', spline_o%ibctype1, ifail)
763  call cdfputvar(ncid, trim(zpre)//'ibctype2', spline_o%ibctype2, ifail)
764  call cdfputvar(ncid, trim(zpre)//'x1', spline_o%x1, ifail)
765  call cdfputvar(ncid, trim(zpre)//'x2', spline_o%x2, ifail)
766  call cdfputvar(ncid, trim(zpre)//'bcval1min', spline_o%bcval1min, ifail)
767  call cdfputvar(ncid, trim(zpre)//'bcval1max', spline_o%bcval1max, ifail)
768  call cdfputvar(ncid, trim(zpre)//'bcval2min', spline_o%bcval2min, ifail)
769  call cdfputvar(ncid, trim(zpre)//'bcval2max', spline_o%bcval2max, ifail)
770  if(ifail/=0) then
771  ier=40
772  return
773  endif
774 
775  if(spline_o%isReady == 1) then
776  if(.not.fullsv) then
777  call cdfputvar(ncid, trim(zpre)//'f', spline_o%fspl(1,:,:), ifail)
778  else
779  call cdfputvar(ncid, trim(zpre)//'x1min', spline_o%x1min, ifail)
780  call cdfputvar(ncid, trim(zpre)//'x1max', spline_o%x1max, ifail)
781  call cdfputvar(ncid, trim(zpre)//'ilin1', spline_o%ilin1, ifail)
782  call cdfputvar(ncid, trim(zpre)//'x2min', spline_o%x2min, ifail)
783  call cdfputvar(ncid, trim(zpre)//'x2max', spline_o%x2max, ifail)
784  call cdfputvar(ncid, trim(zpre)//'ilin2', spline_o%ilin2, ifail)
785  call cdfputvar(ncid, trim(zpre)//'x1pkg', spline_o%x1pkg, ifail)
786  call cdfputvar(ncid, trim(zpre)//'x2pkg', spline_o%x2pkg, ifail)
787  call cdfputvar(ncid, trim(zpre)//'fspl', spline_o%fspl, ifail)
788  endif
789  else
790  ier = 93
791  endif
792 
793  call cdfcls(ncid)
794 
795  if(ifail /=0) ier = 17
796 
797  return
798  end subroutine ezspline_save2_r4
799 
800 
801 
802 
803 
804 !!!
805 !!! 3-D
806 !!!
807 
808 
809 
810  subroutine ezspline_save3_r4(spline_o, filename, ier, &
811  spl_name, fullsave)
812  use ezspline_obj
813  use ezcdf
814  implicit none
815  type(EZspline3_r4) :: spline_o
816  character*(*) :: filename
817  ! ier:
818  ! 17=could not save spline object in file filename
819  integer, intent(out) :: ier
820 
821  ! if SPL_NAME is set, APPEND to file instead of creating new; prepend
822  ! name to all NetCDF data items. This allows one file to contain
823  ! multiple items. (new dmc Mar 2006)
824  character*(*), intent(in), OPTIONAL :: SPL_NAME ! optional spline "name"
825 
826  ! if FULLSAVE is set .TRUE., save derived coefficients along with data
827  ! this saves recalculating the coefficients when file is read (Mar 2006)
828  logical, intent(in), OPTIONAL :: FULLSAVE
829 
830  integer ncid, ifail, in0, in1, in2, in3
831  integer dimlens(3)
832  character(2), parameter :: real4='R4'
833  character(3), parameter :: int='INT'
834 
835  character*32 zpre
836  logical :: fullsv,imodify
837 
838  ier = 0
839  if(spline_o%isReady /= 1) then
840  ier = 94
841  return
842  endif
843 
844  in0 = size(spline_o%fspl,1)
845  in1 = size(spline_o%fspl,2)
846  in2 = size(spline_o%fspl,3)
847  in3 = size(spline_o%fspl,4)
848 
849  if(present(spl_name)) then
850  call ezspline_spl_name_chk(spl_name,ier)
851  if(ier.ne.0) return
852  zpre=spl_name
853  else
854  zpre=' '
855  endif
856 
857  if(present(fullsave)) then
858  fullsv=fullsave
859  else
860  fullsv=.false.
861  endif
862 
863  if(zpre.eq.' ') then
864  call cdfopn(ncid, filename, 'w') ! no error flag??
865  imodify=.false.
866  else
867  call cdfopn(ncid, filename, 'v') ! Append if exists, else create new
868  imodify=.true.
869  endif
870  if(ncid==0) then
871  ier = 39
872  return
873  endif
874 
875  dimlens = (/1, 1, 1/) ! scalar
876  call ezspline_defvar(ncid, trim(zpre)//'n1', dimlens, int, imodify, ier)
877  call ezspline_defvar(ncid, trim(zpre)//'n2', dimlens, int, imodify, ier)
878  call ezspline_defvar(ncid, trim(zpre)//'n3', dimlens, int, imodify, ier)
879  call ezspline_defvar(ncid, trim(zpre)//'klookup1', dimlens, int, imodify, ier)
880  call ezspline_defvar(ncid, trim(zpre)//'klookup2', dimlens, int, imodify, ier)
881  call ezspline_defvar(ncid, trim(zpre)//'klookup3', dimlens, int, imodify, ier)
882  call ezspline_defvar(ncid, trim(zpre)//'isHermite', dimlens, int, imodify, ier)
883  call ezspline_defvar(ncid, trim(zpre)//'isLinear', dimlens, int, imodify, ier)
884  call ezspline_defvar(ncid, trim(zpre)//'isHybrid', dimlens, int, imodify, ier)
885  call ezspline_defvar(ncid, trim(zpre)//'isReady', dimlens, int, imodify, ier)
886  dimlens = (/3, 1, 1/)
887  call ezspline_defvar(ncid, trim(zpre)//'hspline', dimlens, int, imodify, ier)
888  dimlens = (/2, 1, 1/)
889  call ezspline_defvar(ncid, trim(zpre)//'ibctype1', dimlens, int, imodify, ier)
890  call ezspline_defvar(ncid, trim(zpre)//'ibctype2', dimlens, int, imodify, ier)
891  call ezspline_defvar(ncid, trim(zpre)//'ibctype3', dimlens, int, imodify, ier)
892  dimlens = (/spline_o%n1, 1, 1/)
893  call ezspline_defvar(ncid, trim(zpre)//'x1', dimlens, real4, imodify, ier)
894  dimlens = (/spline_o%n2, 1, 1/)
895  call ezspline_defvar(ncid, trim(zpre)//'x2', dimlens, real4, imodify, ier)
896  dimlens = (/spline_o%n3, 1, 1/)
897  call ezspline_defvar(ncid, trim(zpre)//'x3', dimlens, real4, imodify, ier)
898  dimlens = (/in2, in3, 1/)
899  call ezspline_defvar(ncid, trim(zpre)//'bcval1min', dimlens, real4, imodify, ier)
900  call ezspline_defvar(ncid, trim(zpre)//'bcval1max', dimlens, real4, imodify, ier)
901  dimlens = (/in1, in3, 1/)
902  call ezspline_defvar(ncid, trim(zpre)//'bcval2min', dimlens, real4, imodify, ier)
903  call ezspline_defvar(ncid, trim(zpre)//'bcval2max', dimlens, real4, imodify, ier)
904  dimlens = (/in1, in2, 1/)
905  call ezspline_defvar(ncid, trim(zpre)//'bcval3min', dimlens, real4, imodify, ier)
906  call ezspline_defvar(ncid, trim(zpre)//'bcval3max', dimlens, real4, imodify, ier)
907  !
908  ! save only the function values at the nodes
909  !
910  dimlens = (/in1, in2, in3/)
911  if(.not.fullsv) then
912  call ezspline_defvar(ncid, trim(zpre)//'f', dimlens, real4, imodify, ier)
913  else
914  dimlens = (/1, 1, 1/) ! scalar
915  call ezspline_defvar(ncid, trim(zpre)//'x1min', dimlens, real4, imodify, ier)
916  call ezspline_defvar(ncid, trim(zpre)//'x1max', dimlens, real4, imodify, ier)
917  call ezspline_defvar(ncid, trim(zpre)//'ilin1', dimlens, int, imodify, ier)
918  call ezspline_defvar(ncid, trim(zpre)//'x2min', dimlens, real4, imodify, ier)
919  call ezspline_defvar(ncid, trim(zpre)//'x2max', dimlens, real4, imodify, ier)
920  call ezspline_defvar(ncid, trim(zpre)//'ilin2', dimlens, int, imodify, ier)
921  call ezspline_defvar(ncid, trim(zpre)//'x3min', dimlens, real4, imodify, ier)
922  call ezspline_defvar(ncid, trim(zpre)//'x3max', dimlens, real4, imodify, ier)
923  call ezspline_defvar(ncid, trim(zpre)//'ilin3', dimlens, int, imodify, ier)
924  dimlens = (/spline_o%n1, 4, 1/)
925  call ezspline_defvar(ncid, trim(zpre)//'x1pkg', dimlens, real4, imodify, ier)
926  dimlens = (/spline_o%n2, 4, 1/)
927  call ezspline_defvar(ncid, trim(zpre)//'x2pkg', dimlens, real4, imodify, ier)
928  dimlens = (/spline_o%n3, 4, 1/)
929  call ezspline_defvar(ncid, trim(zpre)//'x3pkg', dimlens, real4, imodify, ier)
930  dimlens = (/in0*in1, in2, in3/)
931  call ezspline_defvar(ncid, trim(zpre)//'fspl', dimlens, real4, imodify, ier)
932  endif
933  if(ier.ne.0) return
934 
935  call cdfputvar(ncid, trim(zpre)//'n1', spline_o%n1, ifail)
936  call cdfputvar(ncid, trim(zpre)//'n2', spline_o%n2, ifail)
937  call cdfputvar(ncid, trim(zpre)//'n3', spline_o%n3, ifail)
938  call cdfputvar(ncid, trim(zpre)//'klookup1', spline_o%klookup1, ifail)
939  call cdfputvar(ncid, trim(zpre)//'klookup2', spline_o%klookup2, ifail)
940  call cdfputvar(ncid, trim(zpre)//'klookup3', spline_o%klookup3, ifail)
941  call cdfputvar(ncid, trim(zpre)//'isHermite', spline_o%isHermite, ifail)
942  call cdfputvar(ncid, trim(zpre)//'isLinear', spline_o%isLinear, ifail)
943  call cdfputvar(ncid, trim(zpre)//'isHybrid', spline_o%isHybrid, ifail)
944  call cdfputvar(ncid, trim(zpre)//'isReady', spline_o%isReady, ifail)
945  call cdfputvar(ncid, trim(zpre)//'hspline', spline_o%hspline, ifail)
946  call cdfputvar(ncid, trim(zpre)//'ibctype1', spline_o%ibctype1, ifail)
947  call cdfputvar(ncid, trim(zpre)//'ibctype2', spline_o%ibctype2, ifail)
948  call cdfputvar(ncid, trim(zpre)//'ibctype3', spline_o%ibctype3, ifail)
949  call cdfputvar(ncid, trim(zpre)//'x1', spline_o%x1, ifail)
950  call cdfputvar(ncid, trim(zpre)//'x2', spline_o%x2, ifail)
951  call cdfputvar(ncid, trim(zpre)//'x3', spline_o%x3, ifail)
952  call cdfputvar(ncid, trim(zpre)//'bcval1min', spline_o%bcval1min, ifail)
953  call cdfputvar(ncid, trim(zpre)//'bcval1max', spline_o%bcval1max, ifail)
954  call cdfputvar(ncid, trim(zpre)//'bcval2min', spline_o%bcval2min, ifail)
955  call cdfputvar(ncid, trim(zpre)//'bcval2max', spline_o%bcval2max, ifail)
956  call cdfputvar(ncid, trim(zpre)//'bcval3min', spline_o%bcval3min, ifail)
957  call cdfputvar(ncid, trim(zpre)//'bcval3max', spline_o%bcval3max, ifail)
958  if(ifail/=0) then
959  ier=40
960  return
961  endif
962 
963  if(spline_o%isReady == 1) then
964  if(.not.fullsv) then
965  call cdfputvar(ncid, trim(zpre)//'f', spline_o%fspl(1,:,:,:), ifail)
966  else
967  call cdfputvar(ncid, trim(zpre)//'x1min', spline_o%x1min, ifail)
968  call cdfputvar(ncid, trim(zpre)//'x1max', spline_o%x1max, ifail)
969  call cdfputvar(ncid, trim(zpre)//'ilin1', spline_o%ilin1, ifail)
970  call cdfputvar(ncid, trim(zpre)//'x2min', spline_o%x2min, ifail)
971  call cdfputvar(ncid, trim(zpre)//'x2max', spline_o%x2max, ifail)
972  call cdfputvar(ncid, trim(zpre)//'ilin2', spline_o%ilin2, ifail)
973  call cdfputvar(ncid, trim(zpre)//'x3min', spline_o%x3min, ifail)
974  call cdfputvar(ncid, trim(zpre)//'x3max', spline_o%x3max, ifail)
975  call cdfputvar(ncid, trim(zpre)//'ilin3', spline_o%ilin3, ifail)
976  call cdfputvar(ncid, trim(zpre)//'x1pkg', spline_o%x1pkg, ifail)
977  call cdfputvar(ncid, trim(zpre)//'x2pkg', spline_o%x2pkg, ifail)
978  call cdfputvar(ncid, trim(zpre)//'x3pkg', spline_o%x3pkg, ifail)
979  call ezspline_cdfput3(ncid, trim(zpre)//'fspl', spline_o%fspl, &
980  in0*in1, in2, in3, ifail)
981  endif
982  else
983  ier = 93
984  endif
985 
986  call cdfcls(ncid)
987 
988  if(ifail /=0) ier = 17
989 
990  return
991  end subroutine ezspline_save3_r4
992 
993  subroutine ezspline_spl_name_chk(spl_name,ier)
994  character*(*), intent(in) :: spl_name
995  integer, intent(out) :: ier
996 
997  !----------------
998  integer :: ic,ilen,indx
999  character*63 :: zlegal = &
1000  '0123456789_abcdefghijklmnopqrstuvwxyzABCDEFGHIJKLMNOPQRSTUVWXYZ'
1001  !----------------
1002  ier=0
1003 
1004  ilen = len(trim(spl_name))
1005 
1006  if(ilen.le.0) then
1007  ier=50
1008  else if(ilen.gt.20) then
1009  ier=51
1010  else
1011  indx=index(zlegal,spl_name(1:1))
1012  if(indx.le.10) then
1013  ier=52
1014  else
1015  do ic=2,ilen
1016  indx=index(zlegal,spl_name(ic:ic))
1017  if(indx.le.0) ier=52
1018  enddo
1019  endif
1020  endif
1021 
1022  end subroutine ezspline_spl_name_chk
1023 
1024  subroutine ezspline_defvar(ncid, name, dimlens, ztype, imodify, ier)
1025 
1026  use ezcdf
1027  implicit NONE
1028 
1029  ! inquire about variable with logic specific to ezspline_save
1030 
1031  integer, intent(in) :: ncid ! NetCDF handle
1032  character*(*), intent(in) :: name ! item name to define
1033  integer, intent(in) :: dimlens(3) ! dimensioning information
1034  character*(*), intent(in) :: ztype ! data type e.g. "R8"
1035  logical, intent(in) :: imodify ! T if item may already exist
1036 
1037  integer, intent(inout) :: ier ! status code, set iff a real error
1038  ! is detected...
1039 
1040  !-----------------------------
1041  character*10 :: ztype_loc
1042  integer :: dimlens_loc(3),id1,id2,ii
1043  integer :: ifail,ivarid
1044  !-----------------------------
1045 
1046  if(.not.imodify) then
1047 
1048  ! simply define the quantity...
1049 
1050  call cdfdefvar(ncid, name, dimlens, ztype, ifail)
1051  if(ifail.ne.0) ier=42
1052 
1053  else
1054 
1055  ifail = nf_inq_varid(ncid, name, ivarid)
1056  if(ifail.ne.0) then
1057 
1058  ! assume item does not exist, so, just define it...
1059 
1060  call cdfdefvar(ncid, name, dimlens, ztype, ifail)
1061  if(ifail.ne.0) ier=42
1062 
1063  else
1064 
1065  ! item DOES exist; require dimension & type match...
1066 
1067  call cdfinqvar(ncid, name, dimlens_loc, ztype_loc, ifail)
1068  if(ztype.ne.ztype_loc) then
1069  ier=53
1070 
1071  else
1072  do ii=1,3
1073  id1=max(1,dimlens(ii))
1074  id2=max(1,dimlens_loc(ii))
1075  if(id1.ne.id2) ier=53
1076  enddo
1077  endif
1078 
1079  ! if ier= 0, no cdfDefVar is needed: object already defined
1080  ! and attributes are correct.
1081 
1082  endif
1083  endif
1084 
1085  end subroutine ezspline_defvar
1086 #endif