Scilab Bag Of Tricks: The Scilab-2.5.x IAQ (Infrequently Asked Questions) | ||
---|---|---|
Prev | Chapter 7. Scilab Core | Next |
In the following two sections we shall treat the "anatomy" of native, i.e. low-level Scilab functions. This will confront us with all the gory details of the stack, the low-level API, and the calling conventions. Having the "Guide for Developers", Internals.ps (see also Section 8.2) ready is a good idea. Where the developer guide is at the end of its wits, a study of the source code is appropriate, especially the file SCI/routines/interf/stack1.f is of interest for us.
We start out discussing simple functions in Section 7.3.1. Simple in the sense that they are self-contained and only take non-function parameters as their arguments. In the second part, Section 7.3.2, we shall consider functions that take other functions (either Scilab functions or externals) as arguments. These functionals all rely on correctly set up deisplatch tables, which are treated in Section 7.3.3.
A typical native Scilab function proceeds as follows:
Check the number of input and output parameters.
Get the "pointers" to actual input parameters; supply default values for optional parameters; issue warnings or errors as appropriate if too many or too few parameters are supplied.
Allocate space for temporary variables, "workspace(s)", etc.
It might be necessary to translate the input variables which are in Scilab format into the appropriate format for the worker routine. This is happens for example if the worker routine uses Fortran-77's double complex (or equivalently complex*16) variables. See Section 7.2.2 for details.
Perform the calculations or transformations that really make up the procedure.
As in Step 4, it might be necessary to transform the results, now from the worker routine's format back into Scilab format.
If necessary, allocate space for the return value(s) on the Scilab stack, and copy result(s) to this space.
Now that the general outline is clear, we are ready to dissect a simple function: ortho. It takes exactly one argument a, that is a real or complex m times n matrix. The single output parameter is a matrix of the same shape and type as is the input matrix. The duty of ortho is to transform the columns of the input matrix into orthonormal form; to achieve this we employ the following Lapack functions:
type Complex is record RealPart : Float := 0.0; ImagPart : Float := 0.0; end record; type FloatVector is array (Positive range <>) of Float; type ComplexVector is array (Positive range <>) of Complex; type FloatMatrix is array (1..Lda, Positive range <>) of Float; type ComplexMatrix is array (1..Lda, Positive range <>) of Complex; procedure dgeqrf (M : in Natural; -- number of rows of A N : in Natural; -- number of cols of A A : in out FloatMatrix; -- M-by-N matrix Lda : in Natural; -- leading dimension of A Tau : out FloatVector; -- scalar factors of elementary reflectors Work : out FloatVector; -- workspace Lwork : in Integer; -- size of workspace Work Info : out Integer); -- error indicator procedure dorgqr (M : in Natural; -- number of rows of A N : in Natural; -- number of cols of A K : in Natural; -- number of elementary reflectors A : in out FloatMatrix; -- M-by-N matrix Lda : in Natural; -- leading dimension of A Tau : out FloatVector; -- scalar factors of elementary reflectors Work : out FloatVector; -- workspace Lwork : in Integer; -- size of workspace Work Info : out Integer); -- error indicator procedure zgeqrf (M : in Natural; -- number of rows of A N : in Natural; -- number of cols of A A : in out ComplexMatrix; -- M-by-N matrix Lda : in Natural; -- leading dimension of A Tau : out ComplexVector; -- scalar factors of elementary reflectors Work : out ComplexVector; -- workspace Lwork : in Integer; -- size of workspace Work Info : out Integer); -- error indicator procedure zungqr (M : in Natural; -- number of rows of A N : in Natural; -- number of cols of A K : in Natural; -- number of elementary reflectors A : in out ComplexMatrix; -- M-by-N matrix Lda : in Natural; -- leading dimension of A Tau : out ComplexVector; -- scalar factors of elementary reflectors Work : out ComplexVector; -- workspace Lwork : in Integer; -- size of workspace Work Info : out Integer); -- error indicator procedure dcopy (N : in Natural; -- number of elements to copy X : in FloatVector; -- input vector IncX : in Integer; -- input stride Y : out FloatVector; -- output vector IncY : in Integer); -- output stride
The dgeqrf- and zgeqrf-functions compute a QR-factorization of a real or complex m-by-n matrix a, while the dorgqr-, and zungqr-functions generate an m-by-n real or complex matrix q with orthonormal columns, relying on the QR-factorization of dgeqrf or zgeqrf. Function dcopy copies N elements (of type Float) of the vector X in increments of IncX to the vector Y using increments of IncY on that side. For a detailed description please consult the Lapack User Guide, or the appropriate manual pages. For the mathematics behind the operation, the reader is referred to [Golub:1996].
Example 7-2 is one of the longest examples in the running text, but do not be scared as we will explain line-by-line and variable-by-variable what is there and why.
Example 7-2. Simple native Scilab function
subroutine ortho -- Native functions are parameterless implicit none -- Switch into weeny mode :-) * CONSTANTS integer realtype parameter (realtype = 0) -- See Table 7-2 for type association * LOCAL VARIABLES character*6 fname -- name of the routine as string logical checklhs, checkrhs, cremat, getmat -- Scilab API functions integer topk integer n, m, mattyp integer tausz, worksz, info integer areadr, aimadr, badr, tauadr integer wrkadr, rreadr, rimadr, dumadr * EXTERNAL FUNCTIONS/SUBROUTINES external checklhs, checkrhs, cremat, getmat -- Scilab API functions external error external dcopy, dgeqrf, dorgqr, zgeqrf, zungqr -- LAPACK/BLAS worker subroutines * HEADER include '/site/X11R6/src/scilab/routines/stack.h' -- Scilab API header * TEXT fname = 'ortho' -- Function name (for error messages) topk = top -- top is defined in stack.h rhs = max(0, rhs) if (.not. checkrhs(fname, 1, 1)) returnif (.not. checklhs(fname, 1, 1)) return * fetch input parameters
if (.not. getmat(fname, topk, top - rhs + 1, $ mattyp, m, n, areadr, aimadr)) return if (n * m .eq. 0) return -- Quick return on empty matrix tausz = min(m, n) -- Prescribed by man-page worksz = max(1, n) -- ... same here if (mattyp .eq. realtype) then * real case * allocate temporary variables; all are real
if (.not. cremat(fname, top + 1, realtype, tausz, 1, $ tauadr, dumadr)) return if (.not. cremat(fname, top + 2, realtype, worksz, 1, $ wrkadr, dumadr)) return if (.not. cremat(fname, top + 3, realtype, m, n, $ badr, dumadr)) return * prepare worker routines' input parameters
call dcopy(n * m, stk(areadr), 1, stk(badr), 1) * call worker routines
call dgeqrf(m, n, stk(badr), m, stk(tauadr), $ stk(wrkadr), worksz, info) if (info .ne. 0) then -- Any error is considered fatal buf = fname // ' dgeqrf failed' call error(999) return endif call dorgqr(m, n, tausz, stk(badr), m, stk(tauadr), $ stk(wrkadr), worksz, info) if (info .ne. 0) then -- Any error is considered fatal buf = fname // ' dorgqr failed' call error(999) return endif else * complex case; mattyp != realtype * allocate temporary variables, * use two REAL*8 for one COMPLEX*16
if (.not. cremat(fname, top + 1, realtype, 2 * tausz, 1, $ tauadr, dumadr)) return if (.not. cremat(fname, top + 2, realtype, 2 * worksz, 1, $ wrkadr, dumadr)) return if (.not. cremat(fname, top + 3, realtype, 2 * m, 2 * n, $ badr, dumadr)) return * prepare worker routines' input parameters, joining * two REAL*8 arrays into one COMPLEX*16 array
call dcopy(n * m, stk(areadr), 1, stk(badr), 2) call dcopy(n * m, stk(aimadr), 1, stk(badr + 1), 2) * call worker routines
call zgeqrf(m, n, stk(badr), m, stk(tauadr), $ stk(wrkadr), worksz, info) if (info .ne. 0) then -- Any error is considered fatal buf = fname // ' zgeqrf failed' call error(999) return endif call zungqr(m, n, tausz, stk(badr), m, stk(tauadr), $ stk(wrkadr), worksz, info) if (info .ne. 0) then -- Any error is considered fatal buf = fname // ' zorgqr failed' call error(999) return endif endif * get ready to exit if (lhs .ge. 1) then
if (.not. cremat(fname, top, mattyp, m, n, $ rreadr, rimadr)) return if (mattyp .eq. realtype) then
call dcopy(m * n, stk(badr), 1, stk(rreadr), 1) else call dcopy(m * n, stk(badr), 2, stk(rreadr), 1) call dcopy(m * n, stk(badr + 1), 2, stk(rimadr), 1) endif endif end
Function getmat is called with the second parameter, topk, holding the value of the parameter stack pointer when the control flow entered ortho. This as well as the function name passed in fname is necessary for the cleanup and messaging in case of an error.
The only parameter we use sits on top of the parameter stack for top - rhs + 1 equals top in our case.
On successful return getmat not only sets the data stack addresses areadr, and aimadr of the real and imaginary parts, but also tells us via mattyp whether the matrix is real complex, and via m, and n how large the matrix is.
The following lines directly depend on the sizes passed back form the core interface, calculating the necessary space for two scratch arrays.
We request a purely real storage for each of the three temporary variables, the third parameter being realtype = 0. Therefore, the index for the imaginary part is a dummy index, called dumadr.
Note how the address of the matrices is passed. The idiom is stk(index), where index has been obtained through one of the get*-, or cre*-functions. The memnonic "stk" means data stack.
call dcopy(n * m, stk(areadr), 1, stk(badr), 2) call dcopy(n * m, stk(aimadr), 1, stk(badr + 1), 2)
The first line says: "Copy m times n elements from the first position of the double precision variable stk(areadr) taking each entry (3rd parameter, read stride = 1) to the double complex output variable stk(badr) filling every other entry (5th parameter, write stride = 2)." The second line does almost the same, but starts off writing at the second element stk(badr + 1), therefore filling in the imaginary parts of stk(badr). This corresponds to Step 4.
The situation is a bit more complicated for a complex result, as we have to de-splice the double complex result from Lapack into two double precision matrices. Here are the crucial lines again:
call dcopy(m * n, stk(badr), 2, stk(rreadr), 1) call dcopy(m * n, stk(badr + 1), 2, stk(rimadr), 1)
The first line says: "Copy m times n elements from the first position in the double complex result stk(badr) taking every other entry (3rd parameter, read stride = 2) into the double precision output variable stk(rreadr) filling each entry (5th parameter, write stride = 1)." The second line does almost the same, but starts off at the second element, stk(badr + 1), therfore copying the imaginary parts into stk(rimadr). This way we are merging Step 6, and Step 7 into one.
Func what? What are you talking about? Functionals – what is this? Glad you asked! Functions operate on numbers or variables, which themselves are not functions. The square root function for example is usually applied to numbers (like: sqrt(2)) or more generally to variables (like: sqrt(x) for any real x). Functionals operate on other functions. Prominent examples are differentiation
where f as to fulfill certain continuity requirements at the point x0; integration:
where again f has to fulfill certain (interesting) requirements; and Fourier transformation:
for suitable functions , and integrals
.
The question how to write native Scilab functions that take arbitrary non-function parameters as their arguments has been discussed in the previous section. Now we focus on Scilab functions that take other Scilab functions as their arguments. If the reader does not feel familiar with native Scilab functions, she should reconsider Section 7.3.1.
In a similar manner as in the last section, we introduce an example. The example is taken from our Scilab/Quadpack interface available on the web. Among others it features the integrator dqng for sufficiently smooth functions, which has the following signature:
type SimpleFunctionType is access function(X : in Float) return Float; procedure dqng (Function : in SimpleFunctionType; LowerIntervalEnd : in Float; HigherIntervalEnd : in Float; EpsilonAbsolute : in Float; EpsilonRelative : in Float; Result : out Float; ErrorAbsolute : out Float; NumberOfEvaluations : out Natural; ErrorIndicator : out Natural);
Here comes the complete example.
Example 7-3. Scilab functional
subroutine intsm * * name: intsm.f -- Scilab/F77 interface to QUADPACK's dqng * author: Lydia van Dijk * last rev.: Wed Mar 15 23:49:45 UTC 2000 * scilab ver.: 2.5 * compiler: g77-0.5.25 (Linux 2.3.49) * * Scilab signature: * [res, abs_err, n_eval] = intsm(a, b, f, eps_abs, eps_rel) * * Return Values: * res: value of the integral * abs_err: estimate of the absolute error * n_eval: number of function evaluations * * Arguments (mandatory): * a: lower bound of integral * b: upper bound of integral * f: function to integrate with signature y = f(x), * x, y real scalars * * Arguments (optional): * eps_abs: desired absolute error; default: 0.0 * eps_rel: desired relative error; default: 1e-8 implicit none -- Switch into weeny mode :-) include 'stack.h' common /cintg/ namef -- Name of integrand function external bintg, fintg -- gateways, see Section 7.3.3 external setfintg * LOCAL VARIABLES character*6 namef -- Name of the routine as string character*6 fname -- Name of function to be integrated character*8 errstr logical getexternal, getscalar logical type, cremat integer iadr, sadr, neval, ifail, l, idxf, idxa integer topk, lr, lra, lrb, iipal, dummy double precision epsa, epsr, a, b, val, abserr include 'errnum.inc' -- Error numbers are defined here * STATEMENT FUNCTIONS iadr(l) = l + l - 1 -- Accessor for integers on real*8 stack sadr(l) = l/2 + 1 -- Accessor for real* on integer stack * TEXT fname = 'intsm' -- Name of this function if(rhs .lt. 3 .or. rhs .gt. 5) thencall error(39) return endif topk = top -- Remember stack position * pop optional parameters off the stack if(rhs .eq. 5) then
if (.not. getscalar(fname, topk, top, lr)) return epsr = stk(lr) top = top - 1 else epsr = 1.0d-8 -- Scilab default endif if (rhs .ge. 4) then if (.not. getscalar(fname, topk, top, lr)) return epsa = stk(lr) top = top - 1 else epsa = 0.0d0 -- Scilab default endif * pop mandatory parameters off the stack namef = ' ' -- Fill name-string with 6 spaces type = .false. if (.not. getexternal(fname, topk, top, namef, type, setfintg))
$ return idxf = top -- Remember stack position of function f top = top - 1 if (.not. getscalar(fname, topk, top, lrb)) return b = stk(lrb) top = top - 1 if (.not. getscalar(fname, topk, top, lra)) return a = stk(lra) idxa = top -- Remember stack position of argument a top = topk + 1 -- Reset stack index * call integration routine if (type) then
* compiled external call dqng(fintg, a, b, epsa, epsr, val, abserr, neval, ifail) else
* Scilab macro iipal = iadr(lstk(top)) -- Start building a variable description istk(iipal) = 1 -- ? istk(iipal + 1) = iipal + 2 -- ? istk(iipal + 2) = idxf -- ? istk(iipal + 3) = idxa -- ? lstk(top + 1) = sadr(iipal + 4) -- ? call dqng(bintg, a, b, epsa, epsr, val, abserr, neval, ifail) endif if (ifail .eq. 1) then
buf = fname // ': max. number of steps reached; ' $ // 'integral too difficult for int_sm' call error(emaxdiv) return endif if (ifail .eq. 6) then buf = fname // ': invalid error bounds' call error(ebounds) return endif if (ifail .ne. 0) then * catch all other errors write(errstr, '(I10)') ifail buf = fname // ': unknown error ' // errstr call error(eunknown) return endif * return values #1, and #2 (val, abserr) replace arguments #1, and * #2 (a, b). top = topk - rhs + 1 stk(lra) = val
if (lhs .ge. 2) then top = top + 1 stk(lrb) = abserr endif * return value #3, neval, needs space on the stack if (lhs .ge. 3) then
top = top + 1 if (.not. cremat(fname, top, 0, 1, 1, lrb, dummy)) return stk(lrb) = dble(neval) -- neval is int, stk() is double precision endif end
The case of an external is easy to handle as getexternal has already taken care of initializing the address to be called fintg. A call to setfintg accomplishes this.
FIXME: add text here.
Almost the same holds for the second return value, abserr, though we only can use its slot if there actually is a return variable.
FIXME: write it!
Also called "gateways".
/* * name: quadqack-gw.c -- gateway for all QUADPACK * interface routines * author: Lydia van Dijk * last rev.: Wed Mar 15 02:22:02 UTC 2000 * compiler: pgcc-2.95.2 19991024 (Linux 2.3.49) */ #include <stack-c.h> /* lives in $SCI/routines */ typedef void (*gatef_t) (void); extern void C2F(intals)(void); extern void C2F(intcau)(void); extern void C2F(intexc)(void); extern void C2F(intfou)(void); extern void C2F(intgen)(void); extern void C2F(intgk)(void); extern void C2F(intinf)(void); extern void C2F(intosc)(void); extern void C2F(intsm)(void); static gatef_t gftab[] = { C2F(intals), C2F(intcau), C2F(intexc), C2F(intfou), C2F(intgen), C2F(intgk), C2F(intinf), C2F(intosc), C2F(intsm) }; int C2F(quadpack_gw)(void) { (*gftab[Fin - 1])(); return 0; }
Scilab script part ...
quadpacklibs = ['/site/src/netlib/quadpack/libquadpack-1.0.so', .. '/site/src/netlib/quadpack/intersci/libqpif-1.0.so'] gateway = 'quadpack_gw' // name of the gateway function interfaces = ['intals', 'intcau', 'intexc', 'intfou', .. 'intgen', 'intgk', 'intinf', 'intosc', .. 'intsm'] addinter(quadpacklibs, gateway, interfaces)
The complete example can be found in Section 10.7.