Interface Changes for the naginterfaces Python Package

This page lists backwards-incompatible API changes in naginterfaces up to and including version 30.2.0.0.

Interfaces Changes in version 28.2.0.0

opt - Optimization

  • Solver naginterfaces.library.opt.handle_solve_ipopt() now requires the return vector from its Hessian callback to be updated in place.

    Old definition:

    def hess(x, idf, sigma, lamda, nnzh, inform):
        hx = foo
        return hx, inform
    

    New definition:

    def hess(x, idf, sigma, lamda, hx, inform):
        hx[:] = foo
        return inform
    

    Also in this solver, length argument nnzu has been dropped and the (optional) multipliers u have been added on input for use in a future release. The size of the returned array u will match that supplied for u on input, or be zero if it was not present.

    Old call:

    handle_solve_ipopt(handle, x, nnzu, objfun=...)
    

    New call:

    u = np.empty(nnzu)
    handle_solve_ipopt(handle, x, objfun=..., u=u)
    
  • Solver naginterfaces.library.opt.handle_solve_bxnl() now requires the vector rdx from its lsqgrd callback to be updated in place.

    Old definition:

    def lsqgrd(x, nres, rdx, inform):
        rdx = foo
        return rdx, inform
    

    New definition:

    def lsqgrd(x, nres, rdx, inform):
        rdx[:] = foo
        return inform
    

Interfaces Changes in version 27.3.0.0

opt - Optimization

  • Argument hx in callback lsqhes in naginterfaces.library.opt.handle_solve_bxnl() is now to be modified in place.

    Old definition:

    def lsqhes(x, lamda, inform):
        hx = something
        return hx, inform
    

    New definition:

    def lsqhes(x, lamda, hx, inform):
        hx[:, :] = something
        return inform
    

Interfaces Changes in version 27.1.0.0

anova - Analysis of Variance

  • Values for the maximum number of treatment means to be computed and the maximum number of terms in the analysis of variance table are now computed internally in naginterfaces.library.anova.factorial().

    Old call:

    factorial(y, lfac, nblock, inter, irdf, mterm, maxt)
    

    New call:

    factorial(y, lfac, nblock, inter, irdf)
    

correg - Correlation and Regression Analysis

dot - Inner Products

glopt - Global Optimization

matop - Matrix Operations, Including Inversion

  • Output array w in naginterfaces.library.matop.real_gen_sparse_lu() has been replaced by a single scalar increase which records the growth value of the computed factorization.

    Old call:

    fac = real_gen_sparse_lu(...)
    increase = fac.w[0]
    

    New call:

    fac = real_gen_sparse_lu(...)
    increase = fac.increase
    

mip - Operations Research

  • The interface for naginterfaces.library.mip.iqp_dense() has been improved.

    Old call:

    iqp_dense(a, bl, bu, cvec, h, intvar, istate, xs, strtgy, comm, qphess=...)
    

    New call:

    iqp_dense(bl, bu, h, intvar, xs, strtgy, comm, a=a, cvec=cvec, istate=istate, qphess=...)
    

    in which optionality has been activated for a, cvec and istate and where the returned ax is now None when there were no linear constraints.

mv - Multivariate Methods

omp - OpenMP Utilities

opt - Optimization

  • Some callbacks in naginterfaces.library.opt.handle_solve_bounds_foas(), naginterfaces.library.opt.handle_solve_dfls() and have modified signatures with respect to the inform argument.

    Sample old call:

    objfun = lambda x, _nres: fx
    handle_solve_dfls(...)
    

    New call:

    objfun = lambda x, _nres, inform: (fx, inform)
    handle_solve_dfls(...)
    

    Solver naginterfaces.library.opt.handle_solve_ipopt() also now requires return vectors from its gradient callbacks to be updated in place.

    Old call:

    objfun = lambda x: fx
    objgrd = lambda x, _nnzfd: fdx
    confun = lambda x, _ncnln: gx
    congrd = lambda x, _nnzgd: gdx
    handle_solve_ipopt(...)
    

    New call:

    objfun = lambda x, inform: (fx, inform)
    def objgrd(x, fdx, inform):
        fdx[:] = [...]
        return inform
    confun = lambda x, _ncnln, inform: (gx, inform)
    def congrd(x, gdx, inform):
        gdx[:] = [...]
        return inform
    handle_solve_ipopt(...)
    

    Please refer to one of the specific documents to confirm a routine’s new interface.

  • Argument objfun is now optional in naginterfaces.library.opt.handle_solve_dfno().

    Old call:

    handle_solve_dfno(handle, objfun, x, ...)
    

    New call:

    handle_solve_dfno(handle, x, objfun=objfun, ...)
    
  • Argument maxcal is now optional in naginterfaces.library.opt.lsq_uncon_quasi_deriv_comp().

    Old call:

    lsq_uncon_quasi_deriv_comp(m, selct, lsqfun, maxcal, ..., lsqmon=None, ...)
    

    New call:

    lsq_uncon_quasi_deriv_comp(m, selct, lsqfun, ..., maxcal=maxcal, lsqmon=None, ...)
    
  • The interfaces for naginterfaces.library.opt.lp_solve(), naginterfaces.library.opt.lsq_lincon_solve() and naginterfaces.library.opt.qp_dense_solve() have been improved.

    For naginterfaces.library.opt.lp_solve(); old call:

    lp_solve(nclin, a, bl, bu, cvec, istate, x, comm, io_manager=None)
    

    New call:

    lp_solve(bl, bu, x, comm, a=a, cvec=cvec, istate=istate, io_manager=None)
    

    in which optionality has been activated for a, cvec and istate, nclin is now inferred from a, and where the returned ax is now None when there were no linear constraints.

    For naginterfaces.library.opt.lsq_lincon_solve(); old call:

    lsq_lincon_solve(c, bl, bu, cvec, istate, kx, x, a, b, comm, io_manager=None)
    

    New call:

    lsq_lincon_solve(bl, bu, x, comm, c=c, cvec=cvec, istate=istate, kx=kx, a=a, b=b, io_manager=None)
    

    in which optionality has been activated for c, cvec, istate, kx, a and b, and where the returned b is now None is some situations.

    For naginterfaces.library.opt.qp_dense_solve(); old call:

    qp_dense_solve(nclin, a, bl, bu, cvec, h, istate, x, comm, qphess=None, data=None, io_manager=None)
    

    New call:

    qp_dense_solve(bl, bu, h, x, comm, a=a, cvec=cvec, qphess=None, istate=istate, data=None, io_manager=None)
    

    in which optionality has been activated for a, cvec and istate, nclin is now inferred from a, and where the returned ax is now None when there were no linear constraints.

pde - Partial Differential Equations

quad - Quadrature

rand - Random Number Generators

sparse - Large Scale Linear Systems

sparseig - Large Scale Eigenproblems

  • Argument mb is now required to be zero-sized when the standard eigenvalue problem is being solved in naginterfaces.library.sparseig.real_symm_band_solve(). (This is the default mode of operation.)

    Old call:

    comm = {}
    # The default problem form:
    real_symm_option('Standard', comm)
    ...
    mb = np.empty((2*kl+`ku`+1, n))
    real_symm_band_solve(..., mb, ...)
    

    New call:

    comm = {}
    real_symm_option('Standard', comm)
    ...
    mb = np.empty((0, n))
    real_symm_band_solve(..., mb, ...)
    

tsa - Time Series Analysis

  • Reference vector r in naginterfaces.library.tsa.uni_smooth_exp() has been promoted to a communication structure. Consequently, the structure will now be initialized for you if you pass in an empty dict, and is modified in place.

    Old call:

    mode = 0
    r = np.empty(13 + (0 if itype in (1, 2, 3) else p))
    e = uni_smooth_exp(mode, itype, p, param, y, k, init, nf, r)
    r = e.r
    ...
    

    New call:

    mode = 0
    comm = {}
    e = uni_smooth_exp(mode, itype, p, param, y, k, init, nf, comm)
    ...
    
  • Argument l has been dropped from naginterfaces.library.tsa.multi_kalman_sqrt_invar().

    Old call:

    multi_kalman_sqrt_invar(transf, l, a, b, ...)
    

    New call:

    multi_kalman_sqrt_invar(transf, a, b, ...)
    

Interfaces Changes in version 27.0.0.0

Cross-Chapter Changes

blas - Linear Algebra Support Routines

blgm - Linear Model Specification

  • A new argument what has been added to naginterfaces.library.blgm.lm_submodel() for controlling the labels that are to be produced.

    An old call, in which hform was produced by naginterfaces.library.blgm.lm_​formula(), of the form

    lm_submodel(hform, [...])
    

    must be converted to

    lm_submodel('S', hform, [...])
    

    Additionally, some error exits have been renumbered as follows

    Meaning

    Old code

    New code

    Illegal what

    N/A

    11

    what not valid for the supplied hxdesc

    N/A

    12

    hform not initialized or corrupt

    11

    21

    hform not a G22 handle generated by lm_formula

    12

    22

    A variable name used when creating hform is not present in hxdesc

    13

    23

    Model and design matrix not consistent (w.r.t. mean effect)

    14

    24

    Model and design matrix not consistent (w.r.t. terms)

    15

    25

    Model and design matrix not consistent (w.r.t. mean or main effect)

    16

    26

    Model and design matrix not consistent (w.r.t. contrasts)

    17

    27

    hxdesc not initialized or corrupt

    21

    31

    hxdesc not a G22 handle generated by lm_formula

    22

    32

    Constraint violated for lisx and

    41

    61

    Constraint violated for lplab and

    71

    81

    lenlab too small

    81

    91

    Constraint violated for lvinfo

    91

    101

    lisx, lplab or lvinfo too small

    92

    102

interp - Interpolation

lapackeig - Least Squares and Eigenvalue Problems (LAPACK)

  • Output argument work in naginterfaces.library.lapackeig.dgesvd() has been renamed to workb.

    This affects calls using namedtuple referencing on the output.

    Old call:

    soln = dgesvd([...])
    foo = soln.work
    

    New call:

    soln = dgesvd([...])
    foo = soln.workb
    

matop - Matrix Operations, Including Inversion

mip - Operations Research

  • Argument comm has been dropped from the callback qphess in naginterfaces.library.mip.iqp_dense().

    Old signature:

    qphess(jthcol, h, x, comm, data)
    

    New signature:

    qphess(jthcol, h, x, data)
    

ode - Ordinary Differential Equations

opt - Optimization

pde - Partial Differential Equations

rand - Random Number Generators

  • The population vector ipop is now optional in naginterfaces.library.rand.sample().

    To preserve prior operation

    sample(ipop, m, statecomm)
    

    it is now necessary to use

    sample(m, statecomm, ipop=ipop)
    

    Additional convenience now comes if the input population is , in which case ipop can be omitted and n supplied instead; that is, an existing call of the form

    sample(list(range(1, n+1)), m, statecomm)
    

    should become

    sample(m, statecomm, n=n)
    

    Note that it is not permitted to omit both ipop and n in the new call.

roots - Roots of One or More Transcendental Equations

sparse - Large Scale Linear Systems

stat - Simple Calculations on Statistical Data

  • Argument rcomm in naginterfaces.library.stat.summary_onevar() is now modified in place on the right-hand side of the function call.

    A previous call sequence of the form

    pn = 0
    rcomm = [0.]*20
    while ...:
        summary = summary_onevar(x, ..., pn=pn, rcomm=rcomm)
        pn = summary.pn
        rcomm = summary.comm
    

    should be updated to be

    pn = 0
    rcomm = numpy.empty(20)
    while ...:
        summary = summary_onevar(x, ..., pn=pn, rcomm=rcomm)
        pn = summary.pn
    

sum - Summation of Series

  • The signature for callback f in naginterfaces.library.sum.invlaplace_crump() has been simplified.

    Old signature:

    (fr, fi) = f(pr, pi, data=None)
    

    New signature:

    fp = f(p, data=None)
    

    in which p = complex(pr, pi) and fp = complex(fr, fi).

tsa - Time Series Analysis

wav - Wavelet Transforms

zeros - Zeros of Polynomials

  • Previously-split complex data which appeared as floats in various forms has been coalesced back to complex in the following routines.

    Old call

    zsm, zlg = quadratic_complex(ar, ai, br, bi, cr, ci)
    

    New call

    z_small, z_large = quadratic_complex(a, b, c)
    

    in which a = complex(ar, ai), and so on, and z_small = complex(zsm[0], zsm[1]), etc.

    Old call

    zsm, zlg = quadratic_real(a, b, c)
    

    New call

    z_small, z_large = quadratic_real(a, b, c)
    

    in which z_small = complex(zsm[0], zsm[1]), and so on.

    Old call

    zeror, zeroi, errest = quartic_real(e, a, b, c, d)
    

    New call

    zero, errest = quartic_real(e, a, b, c, d)
    

    in which zero = zeror + zeroi*1j.

    Old call

    zeror, zeroi, errest = quartic_complex(e, a, b, c, d)
    

    New call

    zero, errest = quartic_complex(e, a, b, c, d)
    

    in which zero = zeror + zeroi*1j.

  • Naming in naginterfaces.library.zeros.poly_complex() has been unified with other routines in the submodule.

    This affects calls using named keyword arguments.

    Old call

    poly_complex(ca=ca, scal=scal)
    

    New call

    poly_complex(a=ca, scal=scal)