home *** CD-ROM | disk | FTP | other *** search
/ OS/2 Shareware BBS: 22 gnu / 22-gnu.zip / mesch12a.zip / DOC / tutorial.txt < prev   
Text File  |  1994-01-14  |  47KB  |  1,321 lines

  1.  
  2.  
  3.                MESCHACH VERSION 1.2A
  4.                ---------------------
  5.  
  6.  
  7.                  TUTORIAL
  8.                  ========
  9.  
  10.  
  11.    In this manual the basic data structures are introduced, and some of the
  12. more basic operations are illustrated.  Then some examples of how to use
  13. the data structures and procedures to solve some simple problems are given.
  14. The first example program is a simple 4th order Runge-Kutta solver for
  15. ordinary differential equations.  The second is a general least squares
  16. equation solver for over-determined equations.  The third example
  17. illustrates how to solve a problem involving sparse matrices.  These
  18. examples illustrate the use of matrices, matrix factorisations and solving
  19. systems of linear equations.  The examples described in this manual are
  20. implemented in tutorial.c.
  21.  
  22.    While the description of each aspect of the system is brief and far from
  23. comprehensive, the aim is to show the different aspects of how to set up
  24. programs and routines and how these work in practice, which includes I/O
  25. and error-handling issues.
  26.  
  27.  
  28.  
  29. 1.  THE DATA STRUCTURES AND SOME BASIC OPERATIONS
  30.  
  31.    The three main data structures are those describing vectors, matrices
  32. and permutations.  These have been used to create data structures for
  33. simplex tableaus for linear programming, and used with data structures for
  34. sparse matrices etc.  To use the system reliably, you should always use
  35. pointers to these data structures and use library routines to do all the
  36. necessary initialisation.
  37.  
  38.    In fact, for the operations that involve memory management (creation,
  39. destruction and resizing), it is essential that you use the routines
  40. provided.
  41.  
  42.    For example, to create a matrix A of size 34 , a vector x of dimension
  43. 10, and a permutation p of size 10, use the following code:
  44.  
  45.  
  46.   #include "matrix.h"
  47.   ..............
  48.   main()
  49.   {
  50.      MAT   *A;
  51.      VEC   *x;
  52.      PERM  *p;
  53.      ..........
  54.      A = m_get(3,4);
  55.      x = v_get(10);
  56.      p = px_get(10);
  57.      ..........
  58.   }
  59.  
  60.  
  61.    This initialises these data structures to have the given size.  The
  62. matrix A and the vector x are initially all zero, while p is initially the
  63. identity permutation.
  64.  
  65.    They can be disposed of by calling M_FREE(A), V_FREE(x) and PX_FREE(p)
  66. respectively if you need to re-use the memory for something else.  The
  67. elements of each data structure can be accessed directly using the members
  68. (or fields) of the corresponding structures.  For example the (i,j)
  69. component of A is accessed by A->me[i][j], x_i by x->ve[i] and p_i by
  70. p->pe[i].
  71.  
  72.    Their sizes are also directly accessible: A->m and A->n are the number
  73. of rows and columns of A respectively, x->dim is the dimension of x , and
  74. size of p is p->size.
  75.  
  76.    Note that the indexes are zero relative just as they are in ordinary C,
  77. so that the index i in x->ve[i] can range from 0 to x->dim -1 .  Thus the
  78. total number of entries of a vector is exactly x->dim.
  79.  
  80.    While this alone is sufficient to allow a programmer to do any desired
  81. operation with vectors and matrices it is neither convenient for the
  82. programmer, nor efficient use of the CPU.  A whole library has been
  83. implemented to reduce the burden on the programmer in implementing
  84. algorithms with vectors and matrices.  For instance, to copy a vector from
  85. x to y it is sufficient to write y = v_copy(x,VNULL).  The VNULL is the
  86. NULL vector, and usually tells the routine called to create a vector for
  87. output.
  88.  
  89.    Thus, the v_copy function will create a vector which has the same size
  90. as x and all the components are equal to those of x.  If y has already
  91. been created then you can write y = v_copy(x,y); in general, writing
  92. ``v_copy(x,y);'' is not enough!  If y is NULL, then it is created (to have
  93. the correct size, i.e. the same size as x), and if it is the wrong size,
  94. then it is resized to have the correct size (i.e. same size as x).  Note
  95. that for all the following functions, the output value is returned, even if
  96. you have a non-NULL value as the output argument.  This is the standard
  97. across the entire library.
  98.  
  99.    Addition, subtraction and scalar multiples of vectors can be computed by
  100. calls to library routines: v_add(x,y,out), v_sub(x,y,out), sv_mlt(s,x,out)
  101. where x and y are input vectors (with data type VEC *), out is the output
  102. vector (same data type) and s is a double precision number (data type
  103. double).  There is also a special combination routine, which computes
  104. out=v_1+s,v_2 in a single routine: v_mltadd(v1,v2,s,out).  This is not only
  105. extremely useful, it is also more efficient than using the scalar-vector
  106. multiply and vector addition routines separately.
  107.  
  108.    Inner products can be computed directly: in_prod(x,y) returns the inner
  109. product of x and y.  Note that extended precision evaluation is not
  110. guaranteed.  The standard installation options uses double precision
  111. operations throughout the library.
  112.  
  113.    Equivalent operations can be performed on matrices: m_add(A,B,C) which
  114. returns C=A+B , and sm_mlt(s,A,C) which returns C=sA .  The data types of
  115. A, B and C are all MAT *, while that of s is type double as before.  The
  116. matrix NULL is called MNULL.
  117.  
  118.    Multiplying matrices and vectors can be done by a single function call:
  119. mv_mlt(A,x,out) returns out=A*x while vm_mlt(A,x,out) returns out=A^T*x , or
  120. equivalently, out^T=x^T*A .  Note that there is no distinction between row
  121. and column vectors unlike certain interactive environments such as MATLAB
  122. or MATCALC.
  123.  
  124.    Permutations are also an essential part of the package.  Vectors can be
  125. permuted by using px_vec(p,x,p_x), rows and columns of matrices can be
  126. permuted by using px_rows(p,A,p_A), px_cols(p,A,A_p), and permutations can
  127. be multiplied using px_mlt(p1,p2,p1_p2) and inverted using px_inv(p,p_inv).
  128. The NULL permutation is called PXNULL.
  129.  
  130.    There are also utility routines to initialise or re-initialise these
  131. data structures: v_zero(x), m_zero(A), m_ident(A) (which sets A=I of the
  132. correct size), v_rand(x), m_rand(A) which sets the entries of x and A
  133. respectively to be randomly and uniformly selected between zero and one,
  134. and px_ident(p) which sets p to be an identity permutation.
  135.  
  136.    Input and output are accomplished by library routines v_input(x),
  137. m_input(A), and px_input(p).  If a null object is passed to any of these
  138. input routines, all data will be obtained from the input file, which is
  139. stdin.  If input is taken from a keyboard then the user will be prompted
  140. for all the data items needed; if input is taken from a file, then the
  141. input will have to be of the same format as that produced by the output
  142. routines, which are: v_output(x), m_output(A) and px_output(p).  This
  143. output is both human and machine readable!
  144.  
  145.    If you wish to send the data to a file other than the standard output
  146. device stdout, or receive input from a file or device other than the
  147. standard input device stdin, take the appropriate routine above, use the
  148. ``foutpout'' suffix instead of just ``output'', and add a file pointer as
  149. the first argument.  For example, to send a matrix A to a file called
  150. ``fred'', use the following:
  151.  
  152.  
  153.   #include   "matrix.h"
  154.   .............
  155.   main()
  156.   {
  157.      FILE  *fp;
  158.      MAT   *A;
  159.      .............
  160.      fp = fopen("fred","w");
  161.      m_foutput(fp,A);
  162.      .............
  163.   }
  164.  
  165.  
  166.    These input routines allow for the presence of comments in the data.  A
  167. comment in the input starts with a ``hash'' character ``#'', and continues
  168. to the end of the line.  For example, the following is valid input for a
  169. 3-dimensional vector:
  170.  
  171.   # The initial vector must not be zero
  172.   # x =
  173.   Vector: dim: 3
  174.   -7      0     3
  175.  
  176.  
  177.    For general input/output which conforms to this format, allowing
  178. comments in the input files, use the input() and finput() macros.  These
  179. are used to print out a prompt message if stdin is a terminal (or ``tty''
  180. in Unix jargon), and to skip over any comments if input is from a
  181. non-interactive device.  An example of the usage of these macros is:
  182.  
  183.   input("Input number of steps: ","%d",&steps);
  184.   fp = stdin;
  185.   finput(fp,"Input number of steps: ","%d",&steps);
  186.   fp = fopen("fred","r");
  187.   finput(fp,"Input number of steps: ","%d",&steps);
  188.  
  189. The "%d" is one of the format specifiers which are used in fscanf(); the
  190. last argument is the pointer to the variable (unless the variable is a
  191. string) just as for scanf() and fscanf().  The first two macro calls read
  192. input from stdin, the last from the file fred.  If, in the first two calls,
  193. stdin is a keyboard (a ``tty'' in Unix jargon) then the prompt string
  194.   "Input number of steps: " 
  195. is printed out on the terminal.
  196.  
  197.  
  198.    The second part of the library contains routines for various
  199. factorisation methods.  To use it put
  200.  
  201.   #include   "matrix2.h"
  202.  
  203. at the beginning of your program.  It contains factorisation and solution
  204. routines for LU, Cholesky and QR-factorisation methods, as well as update
  205. routines for Cholesky and QR factorisations.  Supporting these are a number
  206. of Householder transformation and Givens' rotation routines.  Also there is
  207. a routine for generating the Q matrix for a QR-factorisation, if it is
  208. needed explicitly, as it often is.
  209. There are routines for band factorisation and solution for LU and  LDL^T 
  210. factorisations.
  211.  
  212. For using complex numbers, vectors and matrices include
  213.  
  214.   #include   "zmatrix.h"
  215.  
  216. for using the basic routines, and
  217.  
  218.   #include   "zmatrix2.h"
  219.  
  220. for the complex matrix factorisation routines.  The zmatrix2.h file
  221. includes matrix.h and zmatrix.h so you don't need these files included
  222. together.
  223.  
  224. For using the sparse matrix routines in the library you need to put
  225.  
  226.   #include   "sparse.h"
  227.  
  228. or, if you use any sparse factorisation routines,
  229.  
  230.   #include   "sparse2.h"
  231.  
  232. at the beginning of your file.  The routines contained in the library
  233. include routines for creating, destroying, initialising and updating sparse
  234. matrices, and also routines for sparse matrix-dense vector multiplication,
  235. sparse LU factorisation and sparse Cholesky factorisation.
  236.  
  237. For using the iterative routines you need to use
  238.  
  239.   #include   "iter.h"
  240.  
  241. This includes the sparse.h and matrix.h file.
  242. There are also routines for applying iterative methods such as
  243. pre-conditioned conjugate gradient methods to sparse matrices.
  244.  
  245.    And if you use the standard maths library (sin(), cos(), tan(), exp(),
  246. log(), sqrt(), acos() etc.)  don't forget to include the standard
  247. mathematics header:
  248.  
  249.   #include  <math.h>
  250.  
  251. This file is  not  automatically included by any of the Meschach
  252. header files.
  253.  
  254.  
  255.  
  256. 2.  HOW TO MANAGE MEMORY
  257.  
  258.    Unlike many other numerical libraries, Meschach allows you to allocate,
  259. deallocate and resize the vectors, matrices and permutations that you are
  260. using.  To gain maximum benefit from this it is sometimes necessary to
  261. think a little about where memory is allocated and deallocated.  There are
  262. two reasons for this.
  263.  
  264.    Memory allocation, deallocation and resizing takes a significant amount
  265. of time compared with (say) vector operations, so it should not be done too
  266. frequently.  Allocating memory but not deallocating it means that it cannot
  267. be used by any other data structure.  Data structures that are no longer
  268. needed should be explicitly deallocated, or kept as static variables for
  269. later use.  Unlike other interpreted systems (such as Lisp) there is no
  270. implicit ``garbage collection'' of no-longer-used memory.
  271.  
  272.    There are three main strategies that are recommended for deciding how to
  273. allocate, deallocate and resize objects.  These are ``no deallocation''
  274. which is really only useful for demonstration programs, ``allocate and
  275. deallocate'' which minimises overall memory requirements at the expense of
  276. speed, and ``resize on demand'' which is useful for routines that are
  277. called repeatedly.  A new technique for static workspace arrays is to
  278. ``register workspace variables''.
  279.  
  280.  
  281. 2.1  NO DEALLOCATION
  282.  
  283.    This is the strategy of allocating but never deallocating data
  284. structures.  This is only useful for demonstration programs run with small
  285. to medium size data structures.  For example, there could be a line
  286.  
  287.   QR = m_copy(A,MNULL);     /* allocate memory for QR */
  288.  
  289. to allocate the memory, but without the call M_FREE(QR); in it.  This can
  290. be acceptable if QR = m_copy(A,MNULL) is only executed once, and so the
  291. allocated memory never needs to be explicitly deallocated.
  292.  
  293.    This would not be acceptable if QR = m_copy(A,MNULL) occurred inside a
  294. for loop.  If this were so, then memory would be ``lost'' as far as the
  295. program is concerned until there was insufficient space for allocating the
  296. next matrix for QR.  The next subsection shows how to avoid this.
  297.  
  298.  
  299. 2.2  ALLOCATE AND DEALLOCATE
  300.  
  301.    This is the most straightforward way of ensuring that memory is not
  302. lost.  With the example of allocating QR it would work like this:
  303.  
  304.   for ( ... ; ... ; ... )
  305.   {
  306.     QR = m_copy(A,MNULL);    /* allocate memory for QR */
  307.                              /* could have been allocated by m_get() */
  308.     /* use QR */
  309.       ......
  310.       ......
  311.     /* no longer need QR for this cycle */
  312.     M_FREE(QR);             /* deallocate QR so memory can be reused */
  313.   }
  314.  
  315.    The allocate and deallocate statements could also have come at the
  316. beginning and end of a function or procedure, so that when the function
  317. returns, all the memory that the function has allocated has been
  318. deallocated.
  319.  
  320.    This is most suitable for functions or sections of code that are called
  321. repeatedly but involve fairly extensive calculations (at least a
  322. matrix-matrix multiply, or solving a system of equations).
  323.  
  324.  
  325. 2.3  RESIZE ON DEMAND
  326.  
  327.    This technique reduces the time involved in memory allocation for code
  328. that is repeatedly called or used, especially where the same size matrix or
  329. vector is needed.  For example, the vectors v1, v2, etc. in the
  330. Runge-Kutta routine rk4() are allocated according to this strategy:
  331.  
  332.   rk4(...,x,...)
  333.   {
  334.      static VEC *v1=VNULL, *v2=VNULL, *v3=VNULL, *v4=VNULL, *temp=VNULL;
  335.      .......
  336.      v1   = v_resize(v1,x->dim);
  337.      v2   = v_resize(v2,x->dim);
  338.      v3   = v_resize(v3,x->dim);
  339.      v4   = v_resize(v4,x->dim);
  340.      temp = v_resize(temp,x->dim);
  341.      .......
  342.   }
  343.  
  344.    The intention is that the rk4() routine is called repeatedly with the
  345. same size x vector.  It then doesn't make as much sense to allocate v1, v2
  346. etc.  whenever the function is called.  Instead, v_resize() only performs
  347. memory allocation if the memory already allocated to v1, v2 etc. is smaller
  348. than x->dim.
  349.  
  350.    The vectors v1, v2 etc. are declared to be static to ensure that their
  351. values are not lost between function calls.  Variables that are declared
  352. static are set to NULL or zero by default.  So the declaration of v1, v2,
  353. etc., could be
  354.  
  355.   static VEC *v1, *v2, *v3, *v4, *temp;
  356.  
  357.    This strategy of resizing static workspace variables is not so useful if
  358. the object being allocated is extremely large.  The previous ``allocate and
  359. deallocate'' strategy is much more efficient for memory in those
  360. circumstances.  However, the following section shows how to get the best of
  361. both worlds.
  362.  
  363.  
  364. 2.4  REGISTRATION OF WORKSPACE
  365.  
  366.    From version 1.2 onwards, workspace variables can be registered so that
  367. the memory they reference can be freed up on demand.  To do this, the
  368. function containing the static workspace variables has to include calls to
  369. MEM_STAT_REG(var,type) where var is a pointer to a Meschach data type (such
  370. as VEC or MAT).  This call should be placed after the call to the
  371. appropriate resize function.  The type parameter should be a TYPE_... macro
  372. where the ``...'' is the name of a Meschach type such as VEC or MAT.  For
  373. example,
  374.  
  375.   rk4(...,x,...)
  376.   {
  377.      static VEC *v1, *v2, *v3, *v4, *temp;
  378.        .......
  379.      v1   = v_resize(v1,x->dim);
  380.      MEM_STAT_REG(v1,TYPE_VEC);
  381.      v2   = v_resize(v2,x->dim);
  382.      MEM_STAT_REG(v2,TYPE_VEC);
  383.        ......
  384.   }
  385.  
  386. Normally, these registered workspace variables remain allocated.  However,
  387. to implement the ``deallocate on exit'' approach, use the following code:
  388.  
  389.   ......
  390.   mem_stat_mark(1);
  391.   rk4(...,x,...)
  392.   mem_stat_free(1);
  393.   ......
  394.  
  395.    To keep the workspace vectors allocated for the duration of a loop, but
  396. then deallocated, use
  397.  
  398.   ......
  399.   mem_stat_mark(1);
  400.   for (i = 0; i < N; i++ )
  401.     rk4(...,x,...);
  402.   mem_stat_free(1);
  403.   ......
  404.  
  405. The number used in the mem_stat_mark() and mem_stat_free() calls is the
  406. workspace group number.  The call mem_stat_mark(1) designates 1 as the
  407. current workspace group number; the call mem_stat_free(1) deallocates (and
  408. sets to NULL) all static workspace variables registered as belonging to
  409. workspace group 1.
  410.  
  411.  
  412.  
  413. 3.  SIMPLE VECTOR OPERATIONS: AN RK4 ROUTINE
  414.  
  415.    The main purpose of this example is to show how to deal with vectors and
  416. to compute linear combinations.
  417.  
  418.    The problem here is to implement the standard 4th order Runge-Kutta
  419. method for the ODE
  420.  
  421.   x'=f(t,x), x(t_0)=x_0 
  422.  
  423. for x(t_i), i=1,2,3, where t_i=t_0+i*h and h is the step size.
  424.  
  425.    The formulae for the 4th order Runge-Kutta method are:
  426.  
  427.     x_i+1 = x_i+ h/6*(v_1+2*v_2+2*v_3+v_4),
  428. where
  429.     v_1 = f(t_i,x_i)
  430.     v_2 = f(t_i+h, x_i+h*v_1)
  431.     v_3 = f(t_i+h, x_i+h*v_2)
  432.     v_4 = f(t_i+h, x_i+h*v_3)
  433.  
  434. where the v_i are vectors.
  435.  
  436.    The procedure for implementing this method (rk4()) will be passed (a
  437. pointer to) the function f. The implementation of f could, in this system,
  438. create a vector to hold the return value each time it is called.  However,
  439. such a scheme is memory intensive and the calls to the memory allocation
  440. functions could easily dominate the time performed doing numerical
  441. computations.  So, the implementation of f will also be passed an already
  442. allocated vector to be filled in with the appropriate values.
  443.  
  444.    The procedure rk4() will also be passed the current time t, the step
  445. size h, and the current value for x.  The time after the step will be
  446. returned by rk4().
  447.  
  448. The code that does this follows.
  449.  
  450.  
  451.   #include "matrix.h"
  452.  
  453.   /* rk4 - 4th order Runge-Kutta method */
  454.   double rk4(f,t,x,h)
  455.   double t, h;
  456.   VEC    *(*f)(), *x;
  457.   {
  458.      static VEC *v1=VNULL, *v2=VNULL, *v3=VNULL, *v4=VNULL;
  459.      static VEC *temp=VNULL;
  460.  
  461.      /* do not work with NULL initial vector */
  462.      if ( x == VNULL )
  463.         error(E_NULL,"rk4");
  464.  
  465.      /* ensure that v1, ..., v4, temp are of the correct size */
  466.      v1   = v_resize(v1,x->dim);
  467.      v2   = v_resize(v2,x->dim);
  468.      v3   = v_resize(v3,x->dim);
  469.      v4   = v_resize(v4,x->dim);
  470.      temp = v_resize(temp,x->dim);
  471.  
  472.      /* register workspace variables */
  473.      MEM_STAT_REG(v1,TYPE_VEC);
  474.      MEM_STAT_REG(v2,TYPE_VEC);
  475.      MEM_STAT_REG(v3,TYPE_VEC);
  476.      MEM_STAT_REG(v4,TYPE_VEC);
  477.      MEM_STAT_REG(temp,TYPE_VEC);
  478.      /* end of memory allocation */
  479.  
  480.      (*f)(t,x,v1);         /* most compilers allow: f(t,x,v1); */
  481.      v_mltadd(x,v1,0.5*h,temp);    /* temp = x+.5*h*v1 */
  482.      (*f)(t+0.5*h,temp,v2);
  483.      v_mltadd(x,v2,0.5*h,temp);    /* temp = x+.5*h*v2 */
  484.      (*f)(t+0.5*h,temp,v3);
  485.      v_mltadd(x,v3,h,temp);        /* temp = x+h*v3 */
  486.      (*f)(t+h,temp,v4);
  487.  
  488.      /* now add: v1+2*v2+2*v3+v4 */
  489.      v_copy(v1,temp);              /* temp = v1 */
  490.      v_mltadd(temp,v2,2.0,temp);   /* temp = v1+2*v2 */
  491.      v_mltadd(temp,v3,2.0,temp);   /* temp = v1+2*v2+2*v3 */
  492.      v_add(temp,v4,temp);          /* temp = v1+2*v2+2*v3+v4 */
  493.  
  494.      /* adjust x */
  495.      v_mltadd(x,temp,h/6.0,x);     /* x = x+(h/6)*temp */
  496.  
  497.      return t+h;                   /* return the new time */
  498.   }
  499.  
  500.  
  501.    Note that the last parameter of f() is where the output is placed.
  502. Often this can be NULL in which case the appropriate data structure is
  503. allocated and initialised.  Note also that this routine can be used for
  504. problems of arbitrary size, and the dimension of the problem is determined
  505. directly from the data given.  The vectors v_1,...,v_4 are created to have
  506. the correct size in the lines
  507.  
  508.   ....
  509.   v1 = v_resize(v1,x->dim);
  510.   v2 = v_resize(v2,x->dim);
  511.   ....
  512.  
  513.    Here v_resize(v,dim) resizes the VEC structure v to hold a vector of
  514. length dim.  If v is initially NULL, then this creates a new vector of
  515. dimension dim, just as v_get(dim) would do.  For the above piece of code to
  516. work correctly, v1, v2 etc., must be initialised to be NULL vectors.  This
  517. is done by the declaration
  518.  
  519.   static VEC *v1=VNULL, *v2=VNULL, *v3=VNULL, *v4=VNULL;
  520.  
  521. or
  522.  
  523.   static VEC *v1, *v2, *v3, *v4;
  524.  
  525. The operations of vector addition and scalar addition are really the only
  526. vector operations that need to be performed in rk4.  Vector addition is
  527. done by v_add(v1,v2,out), where out=v1+v2, and scalar multiplication by
  528. sv_mlt(scale,v,out), where out=scale*v.
  529.  
  530. These can be combined into a single operation v_mltadd(v1,v2,scale,out),
  531. where out=v1+scale*v2.  As many operations in numerical mathematics involve
  532. accumulating scalar multiples, this is an extremely useful operation, as we
  533. can see above.  For example:
  534.  
  535.   v_mltadd(x,v1,0.5*h,temp);    /* temp = x+0.5*h*v1 */
  536.  
  537.    We also need a number of ``utility'' operations.  For example v_copy(in,
  538. out) copies the vector in to out.  There is also v_zero(v) to zero a vector
  539. v.
  540.  
  541.    Here is an implementation of the function f for simple harmonic motion:
  542.  
  543.   /* f - right-hand side of ODE solver */
  544.   VEC    *f(t,x,out)
  545.   VEC    *x, *out;
  546.   double    t;
  547.   {
  548.     if ( x == VNULL || out == VNULL )
  549.         error(E_NULL,"f");
  550.     if ( x->dim != 2 || out->dim != 2 )
  551.         error(E_SIZES,"f");
  552.  
  553.     out->ve[0] = x->ve[1];
  554.     out->ve[1] = - x->ve[0];
  555.  
  556.     return out;
  557.   }
  558.  
  559.   As can be seen, most of this code is error checking code, which, of
  560. course, makes the routine safer but a little slower.  For a procedure like
  561. f() it is probably not necessary, although then the main program would have
  562. to perform checking to ensure that the vectors involved have the correct
  563. size etc.  The ith component of a vector x is x->ve[i], and indexing is
  564. zero-relative (i.e., the ``first'' component is component 0).  The ODE
  565. described above is for simple harmonic motion:
  566.  
  567.     x_0'=x_1 ,  x_1'=-x_0 , or equivalently,  x_0''+ x_0 = 0 .
  568.  
  569.   Here is the main program:
  570.  
  571.  
  572.   #include <stdio.h>
  573.   #include "matrix.h"
  574.  
  575.   main()
  576.   {
  577.     VEC        *x;
  578.     VEC        *f();
  579.     double     h, t, t_fin;
  580.     double     rk4();
  581.  
  582.     input("Input initial time: ", "%lf", &t);
  583.     input("Input final time: ",  "%lf", &t_fin);
  584.     x = v_get(2);        /* this is the size needed by f() */
  585.     prompter("Input initial state:\n");    x = v_input(VNULL);
  586.     input("Input step size: ", "%lf", &h);
  587.  
  588.     printf("# At time %g, the state is\n",t); 
  589.     v_output(x);
  590.     while ( t < t_fin )
  591.     {
  592.         t = rk4(f,t,x,min(h,t_fin-t));   /* new t is returned */
  593.         printf("# At time %g, the state is\n",t);
  594.         v_output(x);
  595.     t += h;
  596.     }
  597.   }
  598.  
  599.    The initial values are entered as a vector by v_input().  If v_input()
  600. is passed a vector, then this vector will be used to store the input, and
  601. this vector has the size that x had on entry to v_input().  The original
  602. values of x are also used as a prompt on input from a tty.  If a NULL is
  603. passed to v_input() then v_input() will return a vector of whatever size
  604. the user inputs.  So, to ensure that only a two-dimensional vector is used
  605. for the initial conditions (which is what f() is expecting) we use
  606.  
  607.     x = v_get(2);     x = v_input(x);
  608.  
  609.    To compile the program under Unix, if it is in a file tutorial.c:
  610.  
  611.     cc -o tutorial tutorial.c meschach.a
  612.  
  613. or, if you have an ANSI compiler,
  614.  
  615.     cc -DANSI_C -o tutorial tutorial.c meschach.a
  616.  
  617.    Here is a sample session with the above program: 
  618.  
  619.  tutorial
  620.  
  621.   Input initial time: 0
  622.   Input final time: 1
  623.   Input initial state:
  624.   Vector: dim: 2
  625.   entry 0: -1
  626.   entry 1: b
  627.   entry 0: old             -1 new: 1
  628.   entry 1: old              0 new: 0
  629.   Input step size: 0.1
  630.   At time 0, the state is
  631.   Vector: dim: 2
  632.              1              0 
  633.   At time 0.1, the state is
  634.   Vector: dim: 2
  635.     0.995004167  -0.0998333333 
  636.       .................
  637.   At time 1, the state is
  638.   Vector: dim: 2
  639.     0.540302967   -0.841470478 
  640.  
  641.    By way of comparison, the state at t=1 for the true solution is
  642.      x_0(1)=0.5403023058 , x_1(1)=-0.8414709848 .  
  643. The ``b'' that is typed in entering the x vector allows the user to alter
  644. previously entered components. In this case once this is done, the user is
  645. prompted with the old values when entering the new values.  The user can
  646. also type in ``f'' for skipping over the vector's components, which are
  647. then unchanged.  If an incorrectly sized initial value vector x is given,
  648. the error handler comes into action:
  649.  
  650.   Input initial time: 0
  651.   Input final time: 1
  652.   Input initial state:
  653.   Vector: dim: 3
  654.   entry 0: 3
  655.   entry 1: 2
  656.   entry 2: -1
  657.   Input step size: 0.1
  658.   At time 0, the state is
  659.   Vector: dim: 3
  660.              3              2             -1 
  661.  
  662.   "tutorial.c", line 79: sizes of objects don't match in function f()
  663.   Sorry, aborting program
  664.  
  665.    The error handler prints out the error message giving the source code
  666. file and line number as well as the function name where the error was
  667. raised.  The relevant section of f() in file tutorial.c is:
  668.  
  669.   if ( x->dim != 2 || out->dim != 2 )
  670.      error(E_SIZES,"f");               /* line 79 */
  671.  
  672.  
  673.    The standard routines in this system perform error checking of this
  674. type, and also checking for undefined results such as division by zero in
  675. the routines for solving systems of linear equations.  There are also error
  676. messages for incorrectly formatted input and end-of-file conditions.
  677.  
  678.    To round off the discussion of this program, note that we have seen
  679. interactive input of vectors.  If the input file or stream is not a tty
  680. (e.g., a file, a pipeline or a device) then it expects the input to have
  681. the same form as the output for each of the data structures.  Each of the
  682. input routines (v_input(), m_input(), px_input()) skips over ``comments''
  683. in the input data, as do the macros input() and finput().  Anything from a
  684. `#' to the end of the line (or EOF) is considered to be a comment.  For
  685. example, the initial value problem could be set up in a file ivp.dat as:
  686.  
  687.   # Initial time
  688.   0
  689.   # Final time
  690.   1
  691.   # Solution is x(t) = (cos(t),-sin(t))
  692.   # x(0) =
  693.   Vector: dim: 2
  694.   1       0
  695.   # Step size
  696.   0.1
  697.  
  698.    The output of the above program with the above input (from a file) gives
  699. essentially the same output as shown above, except that no prompts are sent
  700. to the screen.
  701.  
  702.  
  703.  
  704. 4.  USING ROUTINES FOR LISTS OF ARGUMENTS
  705.  
  706.    Some of the most common routines have variants that take a variable
  707. number of arguments.  These are the routines .._get_vars(), .._resize_vars()
  708. and .._free_vars().  These correspond to the the basic routines .._get(),
  709. .._resize() and .._free() respectively.  Also there is the
  710. mem_stat_reg_vars() routine which registers a list of static workspace
  711. variables. This corresponds to mem_stat_reg_list() for a single variable.
  712.  
  713.    Here is an example of how to use these functions.  This example also
  714. uses the routine v_linlist() to compute a linear combination of vectors.
  715. Note that the code is much more compact, but don't forget that these
  716. ``..._vars()'' routines usually need the address-of operator ``&'' and NULL
  717. termination of the arguments to work correctly.
  718.  
  719.  
  720.   #include "matrix.h"
  721.  
  722.   /* rk4 - 4th order Runge-Kutta method */
  723.   double rk4(f,t,x,h)
  724.   double t, h;
  725.   VEC    *(*f)(), *x;
  726.   {
  727.     static VEC *v1, *v2, *v3, *v4, *temp;
  728.  
  729.     /* do not work with NULL initial vector */
  730.     if ( x == VNULL )        
  731.     error(E_NULL,"rk4");
  732.  
  733.     /* ensure that v1, ..., v4, temp are of the correct size */
  734.     v_resize_vars(x->dim, &v1, &v2, &v3, &v4, &temp, NULL);
  735.  
  736.     /* register workspace variables */
  737.     mem_stat_reg_vars(0, TYPE_VEC, &v1, &v2, &v3, &v4, &temp, NULL);
  738.     /* end of memory allocation */
  739.  
  740.     (*f)(t,x,v1);             v_mltadd(x,v1,0.5*h,temp);
  741.     (*f)(t+0.5*h,temp,v2);    v_mltadd(x,v2,0.5*h,temp);
  742.     (*f)(t+0.5*h,temp,v3);    v_mltadd(x,v3,h,temp);
  743.     (*f)(t+h,temp,v4);
  744.  
  745.     /* now add: temp = v1+2*v2+2*v3+v4 */
  746.     v_linlist(temp, v1, 1.0, v2, 2.0, v3, 2.0, v4, 1.0, VNULL);
  747.     /* adjust x */
  748.     v_mltadd(x,temp,h/6.0,x);     /* x = x+(h/6)*temp */
  749.  
  750.     return t+h;                   /* return the new time */
  751.   }
  752.  
  753.  
  754.  
  755. 5.  A LEAST SQUARES PROBLEM
  756.  
  757.    Here we need to use matrices and matrix factorisations (in particular, a
  758. QR factorisation) in order to find the best linear least squares solution
  759. to some data.  Thus in order to solve the (approximate) equations
  760.       A*x = b,
  761. where A is an m x n matrix (m > n) we really need to solve the optimisation
  762. problem
  763.       min_x ||Ax-b||^2.  
  764.  
  765.    If we write A=QR where Q is an orthogonal m x m matrix and R is an upper
  766. triangular m x n matrix then (we use 2-norm)
  767.  
  768.     ||A*x-b||^2 = ||R*x-Q^T*b||^2 = || R_1*x - Q_1^T*b||^2 + ||Q_2^T*b||^2
  769.  
  770. where R_1 is an n x n upper triangular matrix.  If A has full rank then R_1
  771. will be an invertible matrix, and the best least squares solution of A*x=b
  772. is x= R_1^{-1}*Q_1^T*b .
  773.  
  774.    These calculations can be be done quite easily as there is a QRfactor()
  775. function available with the system.  QRfactor() is declared to have the
  776. prototype
  777.  
  778.     MAT  *QRfactor(MAT *A, VEC *diag);
  779.  
  780.    The matrix A is overwritten with the factorisation of A ``in compact
  781. form''; that is, while the upper triangular part of A is indeed the R
  782. matrix described above, the Q matrix is stored as a collection of
  783. Householder vectors in the strictly lower triangular part of A and in the
  784. diag vector.  The QRsolve() function knows and uses this compact form and
  785. solves Q*R*x=b with the call QRsolve(A,diag,b,x), which also returns x.
  786.  
  787.    Here is the code to obtain the matrix A, perform the QR factorisation,
  788. obtain the data vector b, solve for x, and determine what the norm of the
  789. errors ( ||Ax-b||_2 ) is.
  790.  
  791.  
  792.   #include "matrix2.h"
  793.  
  794.   main()
  795.   {
  796.     MAT *A, *QR;
  797.     VEC *b, *x, *diag;
  798.  
  799.     /* read in A matrix */
  800.     printf("Input A matrix:");
  801.  
  802.     A = m_input(MNULL);     /* A has whatever size is input */
  803.  
  804.     if ( A->m < A->n )
  805.     {
  806.         printf("Need m >= n to obtain least squares fit");
  807.         exit(0);
  808.     }
  809.     printf("# A =");       m_output(A);
  810.     diag = v_get(A->m);
  811.  
  812.     /* QR is to be the QR factorisation of A */
  813.     QR = m_copy(A,MNULL);
  814.     QRfactor(QR,diag);   
  815.  
  816.     /* read in b vector */
  817.     printf("Input b vector:");
  818.     b = v_get(A->m);
  819.     b = v_input(b);
  820.     printf("# b =");       v_output(b);
  821.  
  822.     /* solve for x */
  823.     x = QRsolve(QR,diag,b,VNULL);
  824.     printf("Vector of best fit parameters is");
  825.     v_output(x);
  826.  
  827.     /* ... and work out norm of errors... */
  828.     printf("||A*x-b|| = %g\n",
  829.     v_norm2(v_sub(mv_mlt(A,x,VNULL),b,VNULL)));
  830.   }
  831.  
  832.    Note that as well as the usual memory allocation functions like m_get(),
  833. the I/O functions like m_input() and m_output(), and the
  834. factorise-and-solve functions QRfactor() and QRsolve(), there are also
  835. functions for matrix-vector multiplication:
  836.      mv_mlt(MAT *A, VEC *x, VEC *out)  
  837. and also vector-matrix multiplication (with the vector on the left):
  838.      vm_mlt(MAT *A, VEC *x, VEC *out), 
  839. with out=x^T A.  There are also functions to perform matrix arithmetic -
  840. matrix addition m_add(), matrix-scalar multiplication sm_mlt(),
  841. matrix-matrix multiplication m_mlt().
  842.  
  843.    Several different sorts of matrix factorisation are supported: LU
  844. factorisation (also known as Gaussian elimination) with partial pivoting,
  845. by LUfactor() and LUsolve().  Other factorisation methods include Cholesky
  846. factorisation CHfactor() and CHsolve(), and QR factorisation with column
  847. pivoting QRCPfactor().
  848.  
  849.    Pivoting involve permutations which have their own PERM data structure.
  850. Permutations can be created by px_get(), read and written by px_input() and
  851. px_output(), multiplied by px_mlt(), inverted by px_inv() and applied to
  852. vectors by px_vec().
  853.  
  854. The above program can be put into a file leastsq.c and compiled under Unix
  855. using
  856.  
  857.     cc -o leastsq leastsq.c meschach.a -lm
  858.  
  859. A sample session using leastsq follows:
  860.  
  861.  
  862.   Input A matrix:
  863.   Matrix: rows cols:5 3
  864.   row 0:
  865.   entry (0,0): 3
  866.   entry (0,1): -1
  867.   entry (0,2): 2
  868.   Continue: 
  869.   row 1:
  870.   entry (1,0): 2
  871.   entry (1,1): -1
  872.   entry (1,2): 1
  873.   Continue: n
  874.   row 1:
  875.   entry (1,0): old              2 new: 2
  876.   entry (1,1): old             -1 new: -1
  877.   entry (1,2): old              1 new: 1.2
  878.   Continue: 
  879.   row 2:
  880.   entry (2,0): old              0 new: 2.5
  881.   ....
  882.   ....             (Data entry)
  883.   ....
  884.   # A =
  885.   Matrix: 5 by 3
  886.   row 0:              3             -1              2 
  887.   row 1:              2             -1            1.2 
  888.   row 2:            2.5              1           -1.5 
  889.   row 3:              3              1              1 
  890.   row 4:             -1              1           -2.2 
  891.   Input b vector:
  892.   entry 0: old              0 new: 5
  893.   entry 1: old              0 new: 3
  894.   entry 2: old              0 new: 2
  895.   entry 3: old              0 new: 4
  896.   entry 4: old              0 new: 6
  897.   # b =
  898.   Vector: dim: 5
  899.            5            3            2            4            6 
  900.   Vector of best fit parameters is
  901.   Vector: dim: 3
  902.      1.47241555   -0.402817858    -1.14411815 
  903.   ||A*x-b|| = 6.78938
  904.  
  905.  
  906.    The Q matrix can be obtained explicitly by the routine makeQ().  The Q
  907. matrix can then be used to obtain an orthogonal basis for the range of A .
  908. An orthogonal basis for the null space of A can be obtained by finding the
  909. QR-factorisation of A^T .
  910.  
  911.  
  912.  
  913. 6.  A SPARSE MATRIX EXAMPLE
  914.  
  915.    To illustrate the sparse matrix routines, consider the problem of
  916. solving Poisson's equation on a square using finite differences, and
  917. incomplete Cholesky factorisation.  The actual equations to solve are
  918.  
  919.     u_{i,j+1} + u_{i,j-1} + u_{i+1,j} + u_{i-1,j} - 4*u_{i,j} =
  920.        h^2*f(x_i,y_j),  for  i,j=1,...,N   
  921.  
  922. where u_{0,j} = u_{i,0} = u_{N+1,j} = u_{i,N+1} = 0 for i,j=1,...,N and h
  923. is the common distance between grid points.
  924.  
  925.    The first task is to set up the matrix describing this system of linear
  926. equations.  The next is to set up the right-hand side.  The third is to
  927. form the incomplete Cholesky factorisation of this matrix, and finally to
  928. use the sparse matrix conjugate gradient routine with the incomplete
  929. Cholesky factorisation as preconditioner.
  930.  
  931.    Setting up the matrix and right-hand side can be done by the following
  932. code:
  933.  
  934.  
  935.   #define N 100
  936.   #define index(i,j) (N*((i)-1)+(j)-1)
  937.   ......
  938.   A = sp_get(N*N,N*N,5);
  939.   b = v_get(N*N);
  940.   h = 1.0/(N+1);      /* for a unit square */
  941.   ......
  942.  
  943.   for ( i = 1; i <= N; i++ )
  944.     for ( j = 1; j <= N; j++ )
  945.     {
  946.         if ( i < N )
  947.             sp_set_val(A,index(i,j),index(i+1,j),-1.0);
  948.         if ( i > 1 )
  949.             sp_set_val(A,index(i,j),index(i-1,j),-1.0);
  950.         if ( j < N )
  951.             sp_set_val(A,index(i,j),index(i,j+1),-1.0);
  952.         if ( j > 1 )
  953.             sp_set_val(A,index(i,j),index(i,j-1),-1.0);
  954.         sp_set_val(A,index(i,j),index(i,j),4.0);
  955.         b->ve[index(i,j)] = -h*h*f(h*i,h*j);
  956.     }
  957.  
  958.    Once the matrix and right-hand side are set up, the next task is to
  959. compute the sparse incomplete Cholesky factorisation of A.  This must be
  960. done in a different matrix, so A must be copied.
  961.  
  962.   LLT = sp_copy(A);
  963.   spICHfactor(LLT);
  964.  
  965. Now when that is done, the remainder is easy:
  966.  
  967.   out = v_get(A->m);
  968.   ......
  969.   iter_spcg(A,LLT,b,1e-6,out,1000,&num_steps);
  970.   printf("Number of iterations = %d\n",num_steps);
  971.   ......
  972.  
  973. and the output can be used in whatever way desired.
  974.  
  975.    For graphical output of the results, the solution vector can be copied
  976. into a square matrix, which is then saved in MATLAB format using m_save(),
  977. and graphical output can be produced by MATLAB.
  978.  
  979.  
  980.  
  981. 7.  HOW DO I ....?
  982.  
  983.    For the convenience of the user, here a number of common tasks that
  984. people need to perform frequently, and how to perform the computations
  985. using Meschach.
  986.  
  987.  
  988. 7.1 .... SOLVE A SYSTEM OF LINEAR EQUATIONS ?
  989.  
  990.    If you wish to solve Ax=b for x given A and b (without destroying A),
  991. then the following code will do this:
  992.  
  993.   VEC   *x, *b;
  994.   MAT    *A, *LU;
  995.   PERM    *pivot;
  996.   ......
  997.   LU = m_get(A->m,A->n);
  998.   LU = m_copy(A,LU);
  999.   pivot = px_get(A->m);
  1000.   LUfactor(LU,pivot);
  1001.   /* set values of b here */
  1002.   x = LUsolve(LU,pivot,b,VNULL);
  1003.  
  1004.  
  1005. 7.2  .... SOLVE A LEAST-SQUARES PROBLEM ?
  1006.  
  1007.    To minimise ||Ax-b||_2^2 = sum_i ((Ax)_i-b_i)^2, the most reliable
  1008. method is based on the QR-factorisation.  The following code performs this
  1009. calculation assuming that A is m x n with m > n :
  1010.  
  1011.   MAT    *A, *QR;
  1012.   VEC    *diag, *b, *x;
  1013.   ......
  1014.   QR = m_get(A->m,A->n);
  1015.   QR = m_copy(A,QR);
  1016.   diag = v_get(A->n);
  1017.   QRfactor(QR,diag);
  1018.   /* set values of b here */
  1019.   x = QRsolve(QR,diag,b,x);
  1020.  
  1021.  
  1022. 7.3  .... FIND ALL THE EIGENVALUES (AND EIGENVECTORS) OF A GENERAL MATRIX ?
  1023.  
  1024.    The best method is based on the Schur decomposition.  For symmetric
  1025. matrices, the eigenvalues and eigenvectors can be computed by a single call
  1026. to symmeig().  For non-symmetric matrices, the situation is more complex
  1027. and the problem of finding eigenvalues and eigenvectors can become quite
  1028. ill-conditioned.  Provided the problem is not too ill-conditioned, the
  1029. following code should give accurate results:
  1030.  
  1031.  
  1032.   /* A is the matrix whose eigenvalues and eigenvectors are sought */
  1033.   MAT    *A, *T, *Q, *X_re, *X_im;
  1034.   VEC    *evals_re, *evals_im;
  1035.   ......
  1036.   Q = m_get(A->m,A->n);
  1037.   T = m_copy(A,MNULL);
  1038.  
  1039.   /* compute Schur form: A = Q*T*Q^T */
  1040.   schur(T,Q);
  1041.  
  1042.   /* extract eigenvalues */
  1043.   evals_re = v_get(A->m);
  1044.   evals_im = v_get(A->m);
  1045.   schur_evals(T,evals_re,evals_im);
  1046.  
  1047.   /* Q not needed for eiegenvalues */
  1048.   X_re = m_get(A->m,A->n);
  1049.   X_im = m_get(A->m,A->n);
  1050.   schur_vecs(T,Q,X_re,X_im);
  1051.   /* k'th eigenvector is k'th column of (X_re + i*X_im) */
  1052.  
  1053.  
  1054.  
  1055. 7.4  .... SOLVE A LARGE, SPARSE, POSITIVE DEFINITE SYSTEM OF EQUATIONS ?
  1056.  
  1057.    An example of a large, sparse, positive definite matrix is the matrix
  1058. obtained from a finite-difference approximation of the Laplacian operator.
  1059. If an explicit representation of such a matrix is available, then the
  1060. following code is suggested as a reasonable way of computing solutions:
  1061.  
  1062.  
  1063.   /* A*x == b is the system to be solved */
  1064.   SPMAT *A, *LLT;
  1065.   VEC    *x, *b;
  1066.   int   num_steps;
  1067.   ......
  1068.   /* set up A and b */
  1069.   ......
  1070.   x = m_get(A->m);
  1071.   LLT = sp_copy(A);
  1072.  
  1073.   /* preconditioning using the incomplete Cholesky factorisation */
  1074.   spICHfactor(LLT);
  1075.  
  1076.   /* now use pre-conditioned conjugate gradients */
  1077.   x = iter_spcg(A,LLT,b,1e-7,x,1000,&num_steps);
  1078.   /* solution computed to give a relative residual of 10^-7 */
  1079.  
  1080.  
  1081.    If explicitly storing such a matrix takes up too much memory, then if
  1082. you can write a routine to perform the calculation of A*x for any given x ,
  1083. the following code may be more suitable (if slower):
  1084.  
  1085.  
  1086.   VEC  *mult_routine(user_def,x,out)
  1087.   void *user_def;
  1088.   VEC  *x, *out;
  1089.   {
  1090.      /* compute out = A*x */
  1091.      ......
  1092.      return out;
  1093.   }
  1094.  
  1095.  
  1096.   main()
  1097.   {
  1098.     ITER *ip;
  1099.     VEC  *x, *b;
  1100.       ......
  1101.     b = v_get(BIG_DIM);     /* right-hand side */
  1102.     x = v_get(BIG_DIM);     /* solution */
  1103.  
  1104.     /* set up b */
  1105.       ......
  1106.     ip = iter_get(b->dim, x->dim);
  1107.     ip->b = v_copy(b,ip->b);
  1108.     ip->info = NULL;        /* if you don't want information
  1109.                                    about solution process */
  1110.     v_zero(ip->x);          /* initial guess is zero */
  1111.     iter_Ax(ip,mult_routine,user_def);
  1112.     iter_cg(ip);
  1113.     printf("# Solution is:\n");   v_output(ip->x);
  1114.       ......
  1115.     ITER_FREE(ip);          /* destroy ip */
  1116.   }
  1117.  
  1118.    The user_def argument is for a pointer to a user-defined structure
  1119. (possibly NULL, if you don't need this) so that you can write a common
  1120. function for handling a large number of different circumstances.
  1121.  
  1122.  
  1123.  
  1124. 8. MORE ADVANCED TOPICS
  1125.  
  1126.    Read this if you are interested in using Meschach library as a base for
  1127. applications. As an example we show how to implement a new type for 3
  1128. dimensional matrices and incorporate this new type into the Meschach
  1129. system. Usually this part of Meschach is transparent to a user.  But a more
  1130. advanced user can take advantage of these routines. We do not describe
  1131. the routines in detail here, but we want to give a rather broad picture of
  1132. what can be done.  By the system we mainly mean the system of delivering
  1133. information on the number of bytes of allocated memory and routines for
  1134. deallocating static variables by mem_stat_... routines.
  1135.  
  1136.    First we introduce a concept of a list of types. By a list of types we
  1137. mean a set of different types with corresponding routines for creating
  1138. these types, destroying and resizing them.  Each type list has a number.
  1139. The list 0 is a list of standard Meschach types such as MAT or VEC. Other
  1140. lists can be defined by a user or a application (based on Meschach). The
  1141. user can attach his/her own list to the system by the routine
  1142. mem_attach_list(). Sometimes it is worth checking if a list number is
  1143. already used by another application. It can be done by
  1144. mem_is_list_attached(ls_num), which returns TRUE if the number ls_num 
  1145. is used. And such a list can be removed from the system by
  1146. mem_free_list(ls_num) if necessary.
  1147.  
  1148.    We describe arguments required by mem_attach_list(). The prototype of
  1149. this function is as follow
  1150.   
  1151.  int mem_attach_list(int ls_num, int ntypes, char *type_names[],
  1152.                  int (*free_funcs[])(), MEM_ARRAY sum[]);
  1153.  
  1154. where the structure MEM_ARRAY has two members: "bytes" of type long and
  1155. "numvar" of type int.  The frst argument is the list number.  Note that you
  1156. cannot overwrite another list.  To do this remove first the old list (by
  1157. mem_free_list()) or choose another number.  The next argument is the number
  1158. of types which are on the list.  This number cannot be changed during
  1159. running a program. The third argument is an array containing the names of
  1160. types (these are character strings).  The fourth one is an array of
  1161. functions deallocating variables of the corresponding type.  And the last
  1162. argument is the local array where information about the number of bytes of
  1163. allocated/deallocated memory (member bytes) and the number of allocated
  1164. variables (member numvar) are gathered. The functions which send
  1165. information to this array are mem_bytes_list() and mem_numvar_list().
  1166.  
  1167.  
  1168. Example:  The routines described here are in the file tutadv.c.
  1169. Firstly we define some macros and a type for 3 dimensional matrices.
  1170.  
  1171. #include "matrix.h"
  1172. #define M3D_LIST    3      /* list number */
  1173. #define TYPE_MAT3D  0      /* the number of a type */
  1174. /* type for 3 dimensional matrices */
  1175. typedef struct {
  1176.     int l,m,n;    /* actual dimensions */
  1177.     int max_l, max_m, max_n;    /* maximal dimensions */
  1178.     Real ***me;    /* pointer to matrix elements */
  1179.                    /* we do not consider segmented memory */
  1180.         Real *base, **me2d;  /* me and me2d are additional pointers 
  1181.                 to base */
  1182. } MAT3D;
  1183.  
  1184.  
  1185. Now we need two routines: one for allocating memory for 3 dimensional
  1186. matrices and the other for deallocating it. It can be useful to have a
  1187. routine for resizing 3 dimensional matrices but we do not use it here.
  1188. Note the use of mem_bytes_list() and mem_numvar_list() to notify the change
  1189. in the number of structures and bytes in use.
  1190.  
  1191. /* function for creating a variable of MAT3D type */
  1192.  
  1193. MAT3D *m3d_get(l,m,n)
  1194. int l,m,n;
  1195. {
  1196.   MAT3D *mat;
  1197.   ....
  1198.   /* alocate memory for structure */
  1199.   if ((mat = NEW(MAT3D)) == (MAT3D *)NULL)
  1200.     error(E_MEM,"m3d_get");
  1201.   else if (mem_info_is_on()) {
  1202.     /* record how many bytes are allocated to structure */
  1203.     mem_bytes_list(TYPE_MAT3D,0,sizeof(MAT3D),M3D_LIST);
  1204.     /* record a new allocated variable */
  1205.     mem_numvar_list(TYPE_MAT3D,1,M3D_LIST);
  1206.   }
  1207.   ....
  1208.   /* allocate memory for 3D array */
  1209.   if ((mat->base = NEW_A(l*m*n,Real)) == (Real *)NULL) 
  1210.     error(E_MEM,"m3d_get");
  1211.   else if (mem_info_is_on())
  1212.     mem_bytes_list(TYPE_MAT3D,0,l*m*n*sizeof(Real),M3D_LIST);
  1213.   ....
  1214.   return mat;
  1215. }
  1216.  
  1217. /* deallocate a variable of type MAT3D */
  1218.  
  1219. int m3d_free(mat)
  1220. MAT3D *mat;
  1221. {
  1222.   /* do not try to deallocate the NULL pointer */
  1223.   if (mat == (MAT3D *)NULL)
  1224.     return -1;
  1225.   ....
  1226.   /* first deallocate base */
  1227.   if (mat->base != (Real *)NULL) {
  1228.     if (mem_info_is_on())
  1229.     /* record how many bytes is deallocated */
  1230.       mem_bytes_list(TYPE_MAT3D,mat->max_l*mat->max_m*mat->max_n*sizeof(Real),
  1231.              0,M3D_LIST);
  1232.     free((char *)mat->base);
  1233.   }
  1234.   ....
  1235.   /* deallocate  MAT3D structure */
  1236.   if (mem_info_is_on()) {
  1237.     mem_bytes_list(TYPE_MAT3D,sizeof(MAT3D),0,M3D_LIST);
  1238.     mem_numvar_list(TYPE_MAT3D,-1,M3D_LIST);
  1239.   }
  1240.   free((char *)mat);
  1241.  
  1242.   ....
  1243.   free((char *)mat);
  1244.  
  1245.   return 0;
  1246. }
  1247.  
  1248.  
  1249. We can now create the arrays necessary for mem_attach_list(). Note that
  1250. m3d_sum can be static if it is in the same file as main(), where
  1251. mem_attach_list is called. Otherwise it must be global.
  1252.  
  1253.  
  1254. char *m3d_names[] = {
  1255.   "MAT3D"
  1256. };
  1257.  
  1258. #define M3D_NUM  (sizeof(m3d_names)/sizeof(*m3d_names))
  1259.  
  1260. int (*m3d_free_funcs[M3D_NUM])() = {
  1261.   m3d_free
  1262. }
  1263.  
  1264. static MEM_ARRAY m3d_sum[M3D_NUM];
  1265.  
  1266.  
  1267. The last thing is to attach the list to the system.
  1268.  
  1269. void main()
  1270. {
  1271.   MAT3D *M;
  1272.   ....
  1273.  
  1274.   mem_info_on(TRUE);    /* switch memory info on */
  1275.   /* attach the new list */
  1276.   mem_attach_list(M3D_LIST,M3D_NUM,m3d_names,m3d_free_funcs,m3d_sum);
  1277.   ....
  1278.   M = m3d_get(3,4,5);
  1279.   ....
  1280.   /* making use of M->me[i][j][k], where i,j,k are non-negative and 
  1281.     i < 3, j < 4, k < 5 */
  1282.   ....
  1283.   mem_info_file(stdout,M3D_LIST);  /* info on the number of allocated 
  1284.                       bytes of memory for types 
  1285.                        on the list M3D_LIST */
  1286.   ....
  1287.   m3d_free(M);  /* if M is not necessary */
  1288.   ....
  1289. }
  1290.  
  1291.  
  1292. We can now use the function mem_info_file() for getting information about
  1293. the number of bytes of allocated memory and number of allocated variables
  1294. of type MAT3D; mem_stat_reg_list() for registering variables of this type
  1295. and mem_stat_mark() and mem_stat_free_list() for deallocating static
  1296. variables of this type.
  1297.  
  1298.  
  1299.  
  1300. In the similar way you can create you own list of errors and attach it to
  1301. the system. See the functions: 
  1302.  
  1303.   int err_list_attach(int list_num, int list_len, char **err_ptr,
  1304.               int warn);  /* for attaching a list of errors */
  1305.  
  1306.   int err_is_list_attached(int list_num);  /* checking if a list 
  1307.                                                     is attached */
  1308.  
  1309.   extern  int err_list_free(int list_num);   /* freeing a list of errors */
  1310.  
  1311. where list_num is the number of the error list, list_len is the number of
  1312. errors on the list, err_ptr is the character string explaining the error
  1313. and warn can be TRUE if this is only a warning (the program continues to
  1314. run) or it can be FALSE if it is an error (the program stops).
  1315.  
  1316. The examples are the standard errors (error list 0) and warnings
  1317. (error list 1) which are in the file err.c
  1318.  
  1319.  
  1320.                 David Stewart and Zbigniew Leyk, 1993
  1321.