~ubuntu-branches/ubuntu/karmic/python-scipy/karmic

« back to all changes in this revision

Viewing changes to Lib/interpolate/fitpack.pyf

  • Committer: Bazaar Package Importer
  • Author(s): Daniel T. Chen (new)
  • Date: 2005-03-16 02:15:29 UTC
  • Revision ID: james.westby@ubuntu.com-20050316021529-xrjlowsejs0cijig
Tags: upstream-0.3.2
ImportĀ upstreamĀ versionĀ 0.3.2

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
!    -*- f90 -*-
 
2
! Author: Pearu Peterson <pearu@cens.ioc.ee>
 
3
!
 
4
python module dfitpack ! in
 
5
 
 
6
  usercode '''
 
7
 
 
8
static double dmax(double* seq,int len) {
 
9
  double val;
 
10
  int i;
 
11
  if (len<1)
 
12
    return -1e308;
 
13
  val = seq[0];
 
14
  for(i=1;i<len;++i)
 
15
    if (seq[i]>val) val = seq[i];
 
16
  return val;
 
17
}
 
18
static double dmin(double* seq,int len) {
 
19
  double val;
 
20
  int i;
 
21
  if (len<1)
 
22
    return 1e308;
 
23
  val = seq[0];
 
24
  for(i=1;i<len;++i)
 
25
    if (seq[i]<val) val = seq[i];
 
26
  return val;
 
27
}
 
28
static double calc_b(double* x,int m,double* tx,int nx) {
 
29
  double val1 = dmin(x,m);
 
30
  double val2 = dmin(tx,nx);
 
31
  if (val2>val1) return val1;
 
32
  val1 = dmax(tx,nx);
 
33
  return val2 - (val1-val2)/nx;
 
34
}
 
35
static double calc_e(double* x,int m,double* tx,int nx) {
 
36
  double val1 = dmax(x,m);
 
37
  double val2 = dmax(tx,nx);
 
38
  if (val2<val1) return val1;
 
39
  val1 = dmin(tx,nx);
 
40
  return val2 + (val2-val1)/nx;
 
41
}
 
42
static int imax(int i1,int i2) {
 
43
  return MAX(i1,i2);
 
44
}
 
45
  '''
 
46
 
 
47
  interface
 
48
 
 
49
     !!!!!!!!!! Univariate spline !!!!!!!!!!!
 
50
 
 
51
     subroutine splev(t,n,c,k,x,y,m,ier)
 
52
       ! y = splev(t,c,k,x)
 
53
       real*8 dimension(n),intent(in) :: t
 
54
       integer intent(hide),depend(t) :: n=len(t)
 
55
       real*8 dimension(n-k-1),depend(n,k),check(len(c)==n-k-1),intent(in) :: c
 
56
       integer :: k
 
57
       real*8 dimension(m),intent(in) :: x
 
58
       real*8 dimension(m),depend(m),intent(out) :: y
 
59
       integer intent(hide),depend(x) :: m=len(x)
 
60
       integer intent(hide) :: ier
 
61
     end subroutine splev
 
62
 
 
63
     subroutine splder(t,n,c,k,nu,x,y,m,wrk,ier)
 
64
       ! dy = splder(t,c,k,x,[nu])
 
65
       real*8 dimension(n) :: t
 
66
       integer depend(t),intent(hide) :: n=len(t)
 
67
       real*8 dimension(n-k-1),depend(n,k),check(len(c)==n-k-1),intent(in) :: c
 
68
       integer :: k
 
69
       integer depend(k),check(0<=nu && nu<=k) :: nu = 1
 
70
       real*8 dimension(m) :: x
 
71
       real*8 dimension(m),depend(m),intent(out) :: y
 
72
       integer depend(x),intent(hide) :: m=len(x)
 
73
       real*8 dimension(n),depend(n),intent(cache,hide) :: wrk
 
74
       integer intent(hide) :: ier
 
75
     end subroutine splder
 
76
 
 
77
     function splint(t,n,c,k,a,b,wrk)
 
78
       ! iy = splint(t,c,k,a,b)
 
79
       real*8 dimension(n),intent(in) :: t
 
80
       integer intent(hide),depend(t) :: n=len(t)
 
81
       real*8 dimension(n),depend(n),check(len(c)==n) :: c
 
82
       integer intent(in) :: k
 
83
       real*8 intent(in) :: a
 
84
       real*8 intent(in) :: b
 
85
       real*8 dimension(n),depend(n),intent(cache,hide) :: wrk
 
86
       real*8 :: splint
 
87
     end function splint
 
88
 
 
89
     subroutine sproot(t,n,c,zero,mest,m,ier)
 
90
       ! zero,m,ier = sproot(t,c[,mest])
 
91
       real*8 dimension(n),intent(in) :: t
 
92
       integer intent(hide),depend(t),check(n>=8) :: n=len(t)
 
93
       real*8 dimension(n),depend(n),check(len(c)==n) :: c
 
94
       real*8 dimension(mest),intent(out),depend(mest) :: zero
 
95
       integer optional,intent(in),depend(n) :: mest=3*(n-7)
 
96
       integer intent(out) :: m
 
97
       integer intent(out) :: ier
 
98
     end subroutine sproot
 
99
 
 
100
     subroutine spalde(t,n,c,k,x,d,ier)
 
101
       ! d,ier = spalde(t,c,k,x)
 
102
 
 
103
       callprotoargument double*,int*,double*,int*,double*,double*,int*
 
104
       callstatement {int k1=k+1; (*f2py_func)(t,&n,c,&k1,&x,d,&ier); }
 
105
 
 
106
       real*8 dimension(n) :: t
 
107
       integer intent(hide),depend(t) :: n=len(t)
 
108
       real*8 dimension(n),depend(n),check(len(c)==n) :: c
 
109
       integer intent(in) :: k
 
110
       real*8 intent(in) :: x
 
111
       real*8 dimension(k+1),intent(out),depend(k) :: d
 
112
       integer intent(out) :: ier
 
113
     end subroutine spalde
 
114
 
 
115
     subroutine fpcurf0(iopt,x,y,w,m,xb,xe,k,s,nest,tol,maxit,k1,k2,n,t,c,fp,fpint,wrk,nrdata,ier)
 
116
       ! x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier = \
 
117
       !   fpcurf0(x,y,k,[w,xb,xe,s,nest])
 
118
 
 
119
       fortranname fpcurf
 
120
       callprotoargument int*,double*,double*,double*,int*,double*,double*,int*,double*,int*,double*,int*,int*,int*,int*,double*,double*,double*,double*,double*,double*,double*,double*,double*,int*,int*
 
121
       callstatement (*f2py_func)(&iopt,x,y,w,&m,&xb,&xe,&k,&s,&nest,&tol,&maxit,&k1,&k2,&n,t,c,&fp,fpint,wrk,wrk+nest,wrk+nest*k2,wrk+nest*2*k2,wrk+nest*3*k2,nrdata,&ier)
 
122
 
 
123
       integer intent(hide) :: iopt = 0
 
124
       real*8 dimension(m),intent(in,out) :: x
 
125
       real*8 dimension(m),depend(m),check(len(y)==m),intent(in,out) :: y
 
126
       real*8 dimension(m),depend(m),check(len(w)==m),intent(in,out) :: w = 1.0
 
127
       integer intent(hide),depend(x),check(m>k),depend(k) :: m=len(x)
 
128
       real*8 intent(in,out),depend(x),check(xb<=x[0]) :: xb = x[0]
 
129
       real*8 intent(in,out),depend(x,m),check(xe>=x[m-1]) :: xe = x[m-1]
 
130
       integer check(1<=k && k<=5),intent(in,out) :: k
 
131
       real*8 check(s>=0.0),depend(m),intent(in,out) :: s = m
 
132
       integer intent(in),depend(m,s,k,k1),check(nest>=2*k1) :: nest = (s==0.0?m+k+1:MAX(m/2,2*k1))
 
133
       real*8 intent(hide) :: tol = 0.001
 
134
       integer intent(hide) :: maxit = 20
 
135
       integer intent(hide),depend(k) :: k1=k+1
 
136
       integer intent(hide),depend(k) :: k2=k+2
 
137
       integer intent(out) :: n
 
138
       real*8 dimension(nest),intent(out),depend(nest) :: t
 
139
       real*8 dimension(nest),depend(nest),intent(out) :: c
 
140
       real*8 intent(out) :: fp
 
141
       real*8 dimension(nest),depend(nest),intent(out,cache)  :: fpint
 
142
       real*8 dimension(nest*3*k2+m*k1),intent(cache,hide),depend(nest,k1,k2,m) :: wrk
 
143
       integer dimension(nest),depend(nest),intent(out,cache) :: nrdata
 
144
       integer intent(out) :: ier
 
145
     end subroutine fpcurf0
 
146
 
 
147
     subroutine fpcurf1(iopt,x,y,w,m,xb,xe,k,s,nest,tol,maxit,k1,k2,n,t,c,fp,fpint,wrk,nrdata,ier)
 
148
       ! x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier = \
 
149
       !   fpcurf1(x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier)
 
150
 
 
151
       fortranname fpcurf
 
152
       callprotoargument int*,double*,double*,double*,int*,double*,double*,int*,double*,int*,double*,int*,int*,int*,int*,double*,double*,double*,double*,double*,double*,double*,double*,double*,int*,int*
 
153
       callstatement (*f2py_func)(&iopt,x,y,w,&m,&xb,&xe,&k,&s,&nest,&tol,&maxit,&k1,&k2,&n,t,c,&fp,fpint,wrk,wrk+nest,wrk+nest*k2,wrk+nest*2*k2,wrk+nest*3*k2,nrdata,&ier)
 
154
 
 
155
       integer intent(hide) :: iopt = 1
 
156
       real*8 dimension(m),intent(in,out,overwrite) :: x
 
157
       real*8 dimension(m),depend(m),check(len(y)==m),intent(in,out,overwrite) :: y
 
158
       real*8 dimension(m),depend(m),check(len(w)==m),intent(in,out,overwrite) :: w
 
159
       integer intent(hide),depend(x),check(m>k),depend(k) :: m=len(x)
 
160
       real*8 intent(in,out) :: xb
 
161
       real*8 intent(in,out) :: xe
 
162
       integer check(1<=k && k<=5),intent(in,out) :: k
 
163
       real*8 check(s>=0.0),intent(in,out) :: s
 
164
       integer intent(hide),depend(t) :: nest = len(t)
 
165
       real*8 intent(hide) :: tol = 0.001
 
166
       integer intent(hide) :: maxit = 20
 
167
       integer intent(hide),depend(k) :: k1=k+1
 
168
       integer intent(hide),depend(k) :: k2=k+2
 
169
       integer intent(in,out) :: n
 
170
       real*8 dimension(nest),intent(in,out,overwrite) :: t
 
171
       real*8 dimension(nest),depend(nest),check(len(c)==nest),intent(in,out,overwrite) :: c
 
172
       real*8 intent(in,out) :: fp
 
173
       real*8 dimension(nest),depend(nest),check(len(fpint)==nest),intent(in,out,cache,overwrite)  :: fpint
 
174
       real*8 dimension(nest*3*k2+m*k1),intent(cache,hide),depend(nest,k1,k2,m) :: wrk
 
175
       integer dimension(nest),depend(nest),check(len(nrdata)==nest),intent(in,out,cache,overwrite) :: nrdata
 
176
       integer intent(in,out) :: ier
 
177
     end subroutine fpcurf1
 
178
 
 
179
     subroutine fpcurfm1(iopt,x,y,w,m,xb,xe,k,s,nest,tol,maxit,k1,k2,n,t,c,fp,fpint,wrk,nrdata,ier)
 
180
       ! x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier = \
 
181
       !   fpcurfm1(x,y,k,t,[w,xb,xe])
 
182
 
 
183
       fortranname fpcurf
 
184
       callprotoargument int*,double*,double*,double*,int*,double*,double*,int*,double*,int*,double*,int*,int*,int*,int*,double*,double*,double*,double*,double*,double*,double*,double*,double*,int*,int*
 
185
       callstatement (*f2py_func)(&iopt,x,y,w,&m,&xb,&xe,&k,&s,&nest,&tol,&maxit,&k1,&k2,&n,t,c,&fp,fpint,wrk,wrk+nest,wrk+nest*k2,wrk+nest*2*k2,wrk+nest*3*k2,nrdata,&ier)
 
186
 
 
187
       integer intent(hide) :: iopt = -1
 
188
       real*8 dimension(m),intent(in,out) :: x
 
189
       real*8 dimension(m),depend(m),check(len(y)==m),intent(in,out) :: y
 
190
       real*8 dimension(m),depend(m),check(len(w)==m),intent(in,out) :: w = 1.0
 
191
       integer intent(hide),depend(x),check(m>k),depend(k) :: m=len(x)
 
192
       real*8 intent(in,out),depend(x),check(xb<=x[0]) :: xb = x[0]
 
193
       real*8 intent(in,out),depend(x,m),check(xe>=x[m-1]) :: xe = x[m-1]
 
194
       integer check(1<=k && k<=5),intent(in,out) :: k
 
195
       real*8 intent(out) :: s = -1
 
196
       integer intent(hide),depend(n) :: nest = n
 
197
       real*8 intent(hide) :: tol = 0.001
 
198
       integer intent(hide) :: maxit = 20
 
199
       integer intent(hide),depend(k) :: k1=k+1
 
200
       integer intent(hide),depend(k) :: k2=k+2
 
201
       integer intent(out),depend(t) :: n = len(t)
 
202
       real*8 dimension(n),intent(in,out,overwrite) :: t
 
203
       real*8 dimension(nest),depend(nest),intent(out) :: c
 
204
       real*8 intent(out) :: fp
 
205
       real*8 dimension(nest),depend(nest),intent(out,cache)  :: fpint
 
206
       real*8 dimension(nest*3*k2+m*k1),intent(cache,hide),depend(nest,k1,k2,m) :: wrk
 
207
       integer dimension(nest),depend(nest),intent(out,cache) :: nrdata
 
208
       integer intent(out) :: ier
 
209
     end subroutine fpcurfm1
 
210
 
 
211
     !!!!!!!!!! Bivariate spline !!!!!!!!!!!
 
212
 
 
213
     subroutine bispev(tx,nx,ty,ny,c,kx,ky,x,mx,y,my,z,wrk,lwrk,iwrk,kwrk,ier)
 
214
       ! z,ier = bispev(tx,ty,c,kx,ky,x,y)
 
215
       real*8 dimension(nx),intent(in) :: tx
 
216
       integer intent(hide),depend(tx) :: nx=len(tx)
 
217
       real*8 dimension(ny),intent(in) :: ty
 
218
       integer intent(hide),depend(ty) :: ny=len(ty)
 
219
       real*8 intent(in),dimension((nx-kx-1)*(ny-ky-1)),depend(nx,ny,kx,ky),&
 
220
            check(len(c)==(nx-kx-1)*(ny-ky-1)):: c
 
221
       integer :: kx
 
222
       integer :: ky
 
223
       real*8 intent(in),dimension(mx) :: x
 
224
       integer intent(hide),depend(x) :: mx=len(x)
 
225
       real*8 intent(in),dimension(my) :: y
 
226
       integer intent(hide),depend(y) :: my=len(y)
 
227
       real*8 dimension(mx,my),depend(mx,my),intent(out,c) :: z
 
228
       real*8 dimension(lwrk),depend(lwrk),intent(hide,cache) :: wrk
 
229
       integer intent(hide),depend(mx,kx,my,ky) :: lwrk=mx*(kx+1)+my*(ky+1)
 
230
       integer dimension(kwrk),depend(kwrk),intent(hide,cache) :: iwrk
 
231
       integer intent(hide),depend(mx,my) :: kwrk=mx+my
 
232
       integer intent(out) :: ier
 
233
     end subroutine bispev
 
234
 
 
235
     subroutine surfit_smth(iopt,m,x,y,z,w,xb,xe,yb,ye,kx,ky,s,nxest,nyest,&
 
236
          nmax,eps,nx,tx,ny,ty,c,fp,wrk1,lwrk1,wrk2,lwrk2,&
 
237
          iwrk,kwrk,ier)
 
238
       !  nx,tx,ny,ty,c,fp,ier = surfit_smth(x,y,z,[w,xb,xe,yb,ye,kx,ky,s,eps,lwrk2])
 
239
 
 
240
       usercode '''
 
241
       static int calc_lwrk1(void) {
 
242
         int u = nxest-kx-1;
 
243
         int v = nyest-ky-1;
 
244
         int km = MAX(kx,ky)+1;
 
245
         int ne = MAX(nxest,nyest);
 
246
         int bx = kx*v+ky+1;
 
247
         int by = ky*u+kx+1;
 
248
         int b1,b2;
 
249
         if (bx<=by) {b1=bx;b2=bx+v-ky;}
 
250
         else {b1=by;b2=by+u-kx;}
 
251
         return u*v*(2+b1+b2)+2*(u+v+km*(m+ne)+ne-kx-ky)+b2+1;
 
252
       }
 
253
       static int calc_lwrk2(void) {
 
254
         int u = nxest-kx-1;
 
255
         int v = nyest-ky-1;
 
256
         int bx = kx*v+ky+1;
 
257
         int by = ky*u+kx+1;
 
258
         int b2 = (bx<=by?bx+v-ky:by+u-kx);
 
259
         return u*v*(b2+1)+b2;
 
260
       }
 
261
       '''
 
262
 
 
263
       fortranname surfit
 
264
       
 
265
       integer intent(hide) :: iopt=0
 
266
       integer intent(hide),depend(x,kx,ky),check(m>=(kx+1)*(ky+1)) &
 
267
            :: m=len(x)
 
268
       real*8 dimension(m) :: x
 
269
       real*8 dimension(m),depend(m),check(len(y)==m) :: y
 
270
       real*8 dimension(m),depend(m),check(len(z)==m) :: z
 
271
       real*8 optional,dimension(m),depend(m),check(len(w)==m) :: w = 1.0
 
272
       real*8 optional,depend(x,m) :: xb=dmin(x,m)
 
273
       real*8 optional,depend(x,m) :: xe=dmax(x,m)
 
274
       real*8 optional,depend(y,m) :: yb=dmin(y,m)
 
275
       real*8 optional,depend(y,m) :: ye=dmax(y,m)
 
276
       integer check(1<=kx && kx<=5) :: kx = 3
 
277
       integer check(1<=ky && ky<=5) :: ky = 3
 
278
       real*8 optional,check(0.0<=s) :: s = m
 
279
       integer optional,depend(kx,m),check(nxest>=2*(kx+1)) &
 
280
            :: nxest = imax(kx+1+sqrt(m/2),2*(kx+1))
 
281
       integer optional,depend(ky,m),check(nyest>=2*(ky+1)) &
 
282
            :: nyest = imax(ky+1+sqrt(m/2),2*(ky+1))
 
283
       integer intent(hide),depend(nxest,nyest) :: nmax=MAX(nxest,nyest)
 
284
       real*8 optional,check(0.0<eps && eps<1.0) :: eps = 1e-16
 
285
       integer intent(out) :: nx
 
286
       real*8 dimension(nmax),intent(out),depend(nmax) :: tx
 
287
       integer intent(out) :: ny
 
288
       real*8 dimension(nmax),intent(out),depend(nmax) :: ty
 
289
       real*8 dimension((nxest-kx-1)*(nyest-ky-1)), &
 
290
            depend(kx,ky,nxest,nyest),intent(out) :: c
 
291
       real*8 intent(out) :: fp
 
292
       real*8 dimension(lwrk1),intent(cache,out),depend(lwrk1) :: wrk1
 
293
       integer intent(hide),depend(m,kx,ky,nxest,nyest) :: lwrk1=calc_lwrk1()
 
294
       real*8 dimension(lwrk2),intent(cache,hide),depend(lwrk2) :: wrk2
 
295
       integer optional,intent(in),depend(kx,ky,nxest,nyest) :: lwrk2=calc_lwrk2()
 
296
       integer dimension(kwrk),depend(kwrk),intent(cache,hide) :: iwrk
 
297
       integer intent(hide),depend(m,nxest,nyest,kx,ky) &
 
298
            :: kwrk=m+(nxest-2*kx-1)*(nyest-2*ky-1)
 
299
       integer intent(out) :: ier
 
300
     end subroutine surfit_smth
 
301
 
 
302
     subroutine surfit_lsq(iopt,m,x,y,z,w,xb,xe,yb,ye,kx,ky,s,nxest,nyest,&
 
303
          nmax,eps,nx,tx,ny,ty,c,fp,wrk1,lwrk1,wrk2,lwrk2,&
 
304
          iwrk,kwrk,ier)
 
305
       ! tx,ty,c,fp,ier = surfit_lsq(x,y,z,tx,ty,[w,xb,xe,yb,ye,kx,ky,eps,lwrk2])
 
306
 
 
307
       usercode '''
 
308
       static int calc_lwrk1(void) {
 
309
         int u = nxest-kx-1;
 
310
         int v = nyest-ky-1;
 
311
         int km = MAX(kx,ky)+1;
 
312
         int ne = MAX(nxest,nyest);
 
313
         int bx = kx*v+ky+1;
 
314
         int by = ky*u+kx+1;
 
315
         int b1,b2;
 
316
         if (bx<=by) {b1=bx;b2=bx+v-ky;}
 
317
         else {b1=by;b2=by+u-kx;}
 
318
         return u*v*(2+b1+b2)+2*(u+v+km*(m+ne)+ne-kx-ky)+b2+1;
 
319
       }
 
320
       static int calc_lwrk2(void) {
 
321
         int u = nxest-kx-1;
 
322
         int v = nyest-ky-1;
 
323
         int bx = kx*v+ky+1;
 
324
         int by = ky*u+kx+1;
 
325
         int b2 = (bx<=by?bx+v-ky:by+u-kx);
 
326
         return u*v*(b2+1)+b2;
 
327
       }
 
328
       '''
 
329
 
 
330
       fortranname surfit
 
331
       
 
332
       integer intent(hide) :: iopt=-1
 
333
       integer intent(hide),depend(x,kx,ky),check(m>=(kx+1)*(ky+1)) &
 
334
            :: m=len(x)
 
335
       real*8 dimension(m) :: x
 
336
       real*8 dimension(m),depend(m),check(len(y)==m) :: y
 
337
       real*8 dimension(m),depend(m),check(len(z)==m) :: z
 
338
       real*8 optional,dimension(m),depend(m),check(len(w)==m) :: w = 1.0
 
339
       real*8 optional,depend(x,tx,m,nx) :: xb=calc_b(x,m,tx,nx)
 
340
       real*8 optional,depend(x,tx,m,nx) :: xe=calc_e(x,m,tx,nx)
 
341
       real*8 optional,depend(y,ty,m,ny) :: yb=calc_b(y,m,ty,ny)
 
342
       real*8 optional,depend(y,ty,m,ny) :: ye=calc_e(y,m,ty,ny)
 
343
       integer check(1<=kx && kx<=5) :: kx = 3
 
344
       integer check(1<=ky && ky<=5) :: ky = 3
 
345
       real*8 intent(hide) :: s = 0.0
 
346
       integer intent(hide),depend(nx) :: nxest = nx
 
347
       integer intent(hide),depend(ny) :: nyest = ny
 
348
       integer intent(hide),depend(nx,ny) :: nmax=MAX(nx,ny)
 
349
       real*8 optional,check(0.0<eps && eps<1.0) :: eps = 1e-16
 
350
       integer intent(hide),depend(tx,kx),check(2*kx+2<=nx) :: nx = len(tx)
 
351
       real*8 dimension(nx),intent(in,out,overwrite) :: tx
 
352
       integer intent(hide),depend(ty,ky),check(2*ky+2<=ny) :: ny = len(ty)
 
353
       real*8 dimension(ny),intent(in,out,overwrite) :: ty
 
354
       real*8 dimension((nx-kx-1)*(ny-ky-1)),depend(kx,ky,nx,ny),intent(out) :: c
 
355
       real*8 intent(out) :: fp
 
356
       real*8 dimension(lwrk1),intent(cache,hide),depend(lwrk1) :: wrk1
 
357
       integer intent(hide),depend(m,kx,ky,nxest,nyest) :: lwrk1=calc_lwrk1()
 
358
       real*8 dimension(lwrk2),intent(cache,hide),depend(lwrk2) :: wrk2
 
359
       integer optional,intent(in),depend(kx,ky,nxest,nyest) :: lwrk2=calc_lwrk2()
 
360
       integer dimension(kwrk),depend(kwrk),intent(cache,hide) :: iwrk
 
361
       integer intent(hide),depend(m,nx,ny,kx,ky) &
 
362
            :: kwrk=m+(nx-2*kx-1)*(ny-2*ky-1)
 
363
       integer intent(out) :: ier
 
364
     end subroutine surfit_lsq
 
365
 
 
366
  end interface
 
367
end python module dfitpack
 
368