SIFEL/GEFEL/matrix.cpp File Reference

#include <stdio.h>
#include <math.h>
#include <string.h>
#include "matrix.h"
#include "eigsol.h"

Go to the source code of this file.

Functions

long addm (const matrix &a, const matrix &b, matrix &c)
 adds 2 matrices C := A+B
long allocim (long n, gfmatrix &mat)
 allocates and initializes identity gfmatrix of dimension n x n
long allocim (long n, matrix &mat)
 allocates and initializes identity matrix of dimension n x n
long allocm (long m, long n, gfmatrix &mat)
 allocates gfmatrix to the dimensions m x n
long allocm (long m, long n, double **&mat)
 allocates double array to the dimensions m x n
long allocm (long m, long n, imatrix &mat)
 allocates imatrix to the dimensions m x n
long allocm (long m, long n, matrix &mat)
 allocates matrix to the dimensions m x n
void bdbj (double *k, const double *b, const double *d, double jac, long m, long n)
 computes matrix product K += B^T . D . B. jac (only upper triangle components are evaluated)
void bdbjac (matrix &d, const matrix &a, const matrix &b, const matrix &c, double jac)
 computes matrix product D += A^T . B . C . jac
long cmulm (double c, matrix &a)
 multiplies matrix by constant A := c.A
long cmulm (const double c, const matrix &a, matrix &b)
 multiplies matrix by constant B := c.A
long condense_matrix (matrix &sm, ivector &cu)
 performs condesation of selected dofs of the given stiffness matrix
long condense_vector (matrix &sm, vector &nf, ivector &cu)
 performs condesation of selected dofs of the given load vector nf according to corresponding stiffness matrix sm
long copym (const matrix &src, matrix &dest)
 copies contents of matrix from the src to dest
long destrm (long **&mat, long m)
 deallocates integer matrix given by 2 dimensional array
long destrm (double **&mat, long m)
 deallocates matrix given by 2 dimensional array
long destrm (matrix &mat)
 deallocates matrix
double det2d (const double *x, const double *y)
 computes determinant of matrix which equals to double area of triangle given by coordinate vectors x and y
double det3d (const double *x, const double *y, const double *z)
 computes determinant of matrix which equals to six-multiple of tetrahedron volume which is given coordinate vectors x,y,z
long detm (const matrix &a, double &det)
 computes determinant of the matrix a, i.e. det(A)
void extractblock (double *a, double *b, long n, long nr, long nc, long fri, long fci)
 extracts submatrix b(nr,nc) from matrix a
void extractm (matrix &a, matrix &b, long fi, long ncomp)
 extracts squared matrix b(ncomp,ncomp) from matrix a
long fillcol (double c, long i, matrix &mat)
 fills i-th column of matrix with given value c
long fillm (double c, matrix &mat)
 fills contents of matrix with given value c
long fillrow (double c, long i, matrix &mat)
 fills i-th row of matrix with given value c
long gause (const matrix &a, matrix &b, vector &r, vector &sol, double zero)
 solves system of linear algebraic equations by the Gauss elimination (A.x = r), resulting matrix is b, resulting vector is stored in r
long gause (const matrix &a, matrix &b, vector &r, double zero)
 solves system of linear algebraic equations by the Gauss elimination (A.x = r), resulting matrix is b, resulting vector is stored in r
long gause (const matrix &a, matrix &b, double zero)
 performs Gauss elimination on the matrix a, result is stored in the matrix b
void gemp (double *a, double *x, double *y, long n, long m, double limit, long pivot)
 solves system of algebraic equations by Gauss method with pivot selection
void glmatrixtransf (matrix &a, const matrix &tmat)
 transforms tensor(matrix) a given in the global coordinate system to the local one (A_l := T^T . A_g . T)
void globloc (const double *gv, double *lv, const long *cn, long n)
 localizes components of local(element) vector to the global(problem) vector
void glvectortransf (const vector &g, vector &l, const matrix &tmat)
 transforms vector g given in the global coordinate system to vector l in local coordinate system (l = T^T . g)
void glvectortransfblock (vector &v, const matrix &t)
 transforms block vector v given in the global coordinate system to the local one (v_l := T^T . v_g) (number of blocks is given by dimension of T)
long identm (matrix &mat)
 sets contents of mat to identity matrix
long invm (const matrix &a, matrix &b, double zero)
 inverts matrix a, i.e. computes A^{-1}
void lgmatrixtransf (matrix &a, const matrix &tmat)
 transforms tensor(matrix) a given in the local coordinate system to the global one (A_g := T . A_l . T^T)
void lgmatrixtransfblock (matrix &a, const matrix &t)
 transforms block matrix a given in the local coordinate system to the global one (A_g := T . A_l . T^T) (number of blocks is given by dimension of T)
void lgvectortransf (vector &g, const vector &l, const matrix &tmat)
 transforms vector l given in the local coordinate system to vector g in global coordinate system (g = T . l)
void lgvectortransfblock (vector &v, const matrix &t)
 transforms block vector v given in the local coordinate system to the global one (v_g := T^T . v_l) (number of blocks is given by dimension of T)
void locglob (double *gv, const double *lv, const long *cn, long n)
 localizes components of local(element) vector to the global(problem) vector
void locvecmat (double *mat, const double *vect, const long *cn, long ci, long m, long n)
void lu_full (double *a, double *x, double *y, long n, double zero, long tc)
 solves system of algebraic equations by LU decomposition
void mat_localize (matrix &gm, matrix &lm, long *rcn, long *ccn)
 localizes components of local matrix lm to the global one gm
void matassem_lsm (double *lsm, vector &natcoord)
 assembles matrix of least square problem for extrapolation of nodal values from integration point values
long mixv (const matrix &a, const vector &u, double &vi, long i)
 computes i-th component of the vector v given by product of matrix a and vector u v := A.u
long mswapc (matrix &mat, long i, long j)
 swaps matrix columns i and j
long mswapr (matrix &mat, long i, long j)
 swaps matrix rows i and j
void mtxm (const double *a, const double *b, double *c, long l, long m, long n)
 multiplies two matrices C := A^T.B given by arrays
long mtxm (const matrix &a, const matrix &b, matrix &c)
 multiplies two matrices C := A^T.B
void mtxmccr (const double *a, const double *b, double *c, long l, long m, long n)
 multiplies transposed matrix a from right by matrix b, C := A^T . B, matrices are stored by columns
void mtxv (const double *a, const double *b, double *c, long m, long n)
 multiplies transposed matrix a given by array from right by vector u given by array v := A^T.u
long mtxv (const matrix &a, const vector &u, vector &v)
 multiplies transposed matrix a from right by vector u v := A^T.u
void mtxvc (const double *a, const double *b, double *c, long m, long n)
void mxm (const double *a, const double *b, double *c, long l, long m, long n)
 multiplies two matrices C := A.B given by arrays
long mxm (const matrix &a, const matrix &b, matrix &c)
 multiplies two matrices C := A.B
void mxmt (const double *a, const double *b, double *c, long l, long m, long n)
 multiplies two matrices C := A.B^T given by arrays
long mxmt (const matrix &a, const matrix &b, matrix &c)
 multiplies two matrices C := A.B^T
void mxv (const double *a, const double *b, double *c, long m, long n)
 multiplies matrix a given by array from right by vector u given by array v := A.u
long mxv (const matrix &a, const vector &u, vector &v)
 multiplies matrix a from right by vector u v := A.u
void mxvc (const double *a, const double *b, double *c, long m, long n)
 multiplies matrix a stored by columns from the right by vector b c := A.b
void nnj (double *m, const double *n, double jac, long mm, long nn)
 computes matrix product M += N^T . N . jac
void nnjac (matrix &d, const matrix &a, const matrix &b, double jac)
 computes matrix product D += A^T . B . jac
void princ_val (matrix &v, vector &pval, matrix &pvect, long ni, double err, double zero, long n, long normalize)
 computes eigenvalues(principal values) of the given matrix v by Jacobi method
long printm (FILE *out, const gfmatrix &a)
 prints the gfmatrix to the file
long printm (FILE *out, const matrix &a)
 prints the matrix to the file
long printm (const matrix &a, FILE *out, int prec, int width)
 prints the matrix to the file with the given width and precision
long readm (XFILE *in, gfmatrix &a)
 reads gfmatrix from the file in
long readm (XFILE *in, matrix &a)
 reads matrix from the file in
long reallocm (long m, long n, imatrix &mat)
 reallocates imatrix to the dimensions m x n
long reallocm (long m, long n, matrix &mat)
 reallocates matrix to the dimensions m x n
void rhsassem_lsm (double *rhs, vector &natcoord, vector &values)
 assembles right hand side vector of least square problem for extrapolation of nodal values from integration point values
void solve_lsm (double *lsm, double *lhs, double *rhs, double zero, long n, long m)
 solves system of equation given by least square problem for extrapolation of nodal values from integration point values
long subm (const matrix &a, const matrix &b, matrix &c)
 substracts 2 matrices C := A-B
long tensprd (const vector &a, const vector &b, matrix &c)
long tranm (matrix &a)
 transposes matrix A := A^T
long tranm (const matrix &a, matrix &at)
 transposes matrix A_t := A^T
long vxm (const vector &u, const matrix &a, vector &v)
 multiplies matrix a from left by row vector u v := u.A
long vxmt (const vector &u, const matrix &a, vector &v)
 multiplies transposed matrix a from left by row vector u v := u.A^T
long vxmxv (const double *v, const double *m, long dim, double &answer)
 multiplies row vector v given by array with symetric matrix m given by array and column vector v, i.e. answer := (v^T).A.v
long vxmxv (const vector &v, const matrix &m, double &answer)
 multiplies row vector v with symetric matrix m and column vector v, i.e. answer := (v^T).A.v
long vxv (const vector &u, const vector &v, matrix &a)
 multiplies column vector u with row vector, i.e. it performs tensor product, A := u * v

Function Documentation

long addm ( const matrix a,
const matrix b,
matrix c 
)

adds 2 matrices C := A+B

The function adds matrix given by a to matrix given by b, the result is stored in c.

Parameters:
a is the structure of the first added matrix
b is the structure of the second added matrix
c is the structure of the result matrix

Requests : a, b and c have to be same dimensions, c has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes.
1 : in case incompatibility sizes of a, b and c matrices.

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 729 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

long allocim ( long  n,
gfmatrix mat 
)

allocates and initializes identity gfmatrix of dimension n x n

The function allocates memory from heap to the mat's member array a and it is initialized to be identity gfmatrix, i.e. all matrix component general functions are stat type with 0 value except of diagonal elements whose value is set to 1 .

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated gfmatrix

created 27.11.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 416 of file matrix.cpp.

References gfmatrix::a, gfunct::f, gfmatrix::m, gfmatrix::n, and print_err().

long allocim ( long  n,
matrix mat 
)

allocates and initializes identity matrix of dimension n x n

The function allocates memory from heap to the mat's member array a and it is initialized to be identity matrix.

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated matrix

created 27.11.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 383 of file matrix.cpp.

References matrix::a, matrix::m, memset(), matrix::n, and print_err().

Referenced by outdriverm::read().

long allocm ( long  m,
long  n,
gfmatrix mat 
)

allocates gfmatrix to the dimensions m x n

The function allocates memory from heap to the mat's member array a.

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated gfmatrix

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 358 of file matrix.cpp.

References gfmatrix::a, gfmatrix::m, gfmatrix::n, and print_err().

long allocm ( long  m,
long  n,
double **&  mat 
)

allocates double array to the dimensions m x n

The function allocates memory from heap to the mat's member array a.

Parameters:
m is the number of rows
n is the number of columns
mat is the allocated double array

created 20.11.2003, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 332 of file matrix.cpp.

References print_err().

long allocm ( long  m,
long  n,
imatrix mat 
)

allocates imatrix to the dimensions m x n

The function allocates memory from heap to the mat's member array a.

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated imatrix

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 304 of file matrix.cpp.

References imatrix::a, imatrix::m, memset(), imatrix::n, and print_err().

long allocm ( long  m,
long  n,
matrix mat 
)

allocates matrix to the dimensions m x n

The function allocates memory from heap to the mat's member array a.

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated matrix

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 275 of file matrix.cpp.

References matrix::a, matrix::m, memset(), matrix::n, and print_err().

void bdbj ( double *  k,
const double *  b,
const double *  d,
double  jac,
long  m,
long  n 
)

computes matrix product K += B^T . D . B. jac (only upper triangle components are evaluated)

The function computes matrix multiplication in the form B^T . D . B . jac and the result is added to the matrix k. It is supposed that the resulting matrix k is symmetric and therefore only upper triangle components are evaluated while the remaining are assigned only.

Parameters:
k - array of resulting matrix K(n,n) = B^T . D . B . jac
b - array of matrix B(m,n)
d - array of matrix D(m,m)
jac - jacobian (constant)
m - number of rows of matrix B, number of rows/columns of matrix D
n - number of columns of matrix B

Created by JK 12.7.1996

Definition at line 2524 of file matrix.cpp.

References ll.

void bdbjac ( matrix d,
const matrix a,
const matrix b,
const matrix c,
double  jac 
)

computes matrix product D += A^T . B . C . jac

The function computes matrix multiplication in the form A^T . B . C . jac and the result is added to the matrix d.

Parameters:
d - resulting matrix D(n,n) += A^T . B . C . jac
a - matrix A(m,n)
b - matrix B(m,m)
b - matrix C(m,n)
jac - jacobian (constant)

Created by JK 12.7.1996

10.5.2002

Definition at line 2587 of file matrix.cpp.

References addm(), allocm(), cmulm(), destrm(), matrix::m, mtxm(), mxm(), and matrix::n.

long cmulm ( double  c,
matrix a 
)

multiplies matrix by constant A := c.A

The function multiplies matrix given by a with constant c, the result is stored in a.

Parameters:
c is the real number
a is the structure of the matrix, from which is multiplied
Returns:
always zero

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1145 of file matrix.cpp.

References matrix::m, and matrix::n.

long cmulm ( const double  c,
const matrix a,
matrix b 
)

multiplies matrix by constant B := c.A

The function multiplies matrix given by a by constant c, the result is stored in b

Parameters:
c is the real number
a is the structure of the matrix, from which is multiplied
b is the structure of the result matrix

Requests : a and b have to be same dimensions, b has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a and b matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1114 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

long condense_matrix ( matrix sm,
ivector cu 
)

performs condesation of selected dofs of the given stiffness matrix

The function performs a condensation of selected dofs in the stiffness matrix sm.

Parameters:
sm - stiffness matrix
cu - vector with number one on the place of condensed dofs

Created 23.5.2012, Jan Zitny, zitny.ja@gmail.com

Definition at line 2409 of file matrix.cpp.

References allocv(), copym(), and matrix::m.

long condense_vector ( matrix sm,
vector nf,
ivector cu 
)

performs condesation of selected dofs of the given load vector nf according to corresponding stiffness matrix sm

The function performs a condensation of selected dofs of the load vector according to corresponding stiffness matrix.

Parameters:
sm - stiffness matrix
f - load vector
cu - is a vector with number one on the position of condensed dofs

Created 23.5.2012, Jan Zitny, zitny.ja@gmail.com

Definition at line 2466 of file matrix.cpp.

References allocv(), and matrix::m.

long copym ( const matrix src,
matrix dest 
)

copies contents of matrix from the src to dest

The function copies matrix given by src to dest.

Parameters:
src is the structure of source matrix to copy
dest is the structure of destination matrix to which will be copied contents of src

Requests : dest has to be setuped dimensions and allocated memory array for elements which is enough large to hold all contents of the matrix src.

Return values:
0 : on succes
1 : in case incompatibility sizes of the src and dest matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 528 of file matrix.cpp.

References matrix::a, matrix::m, matrix::n, and print_err().

long destrm ( long **&  mat,
long  m 
)

deallocates integer matrix given by 2 dimensional array

The function deallocates memory occupied by long array mat.

Parameters:
m - the number of rows
mat - array
Returns:
always zero

created 5.5.2003, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 700 of file matrix.cpp.

long destrm ( double **&  mat,
long  m 
)

deallocates matrix given by 2 dimensional array

The function deallocates memory occupied by double array mat.

Parameters:
m - the number of rows
mat - array
Returns:
always zero

created 5.5.2003, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 678 of file matrix.cpp.

long destrm ( matrix mat  ) 

deallocates matrix

The function deallocates memory occupied by mat's member array a.

Parameters:
mat is the structure for the allocated matrix
Returns:
always zero

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 653 of file matrix.cpp.

References matrix::a, matrix::m, and matrix::n.

double det2d ( const double *  x,
const double *  y 
)

computes determinant of matrix which equals to double area of triangle given by coordinate vectors x and y

The function computes determinant of the matrix

| 1 x1 y1 | | 1 x2 y2 | | 1 x3 y3 |

which is equal to double area of the triangle given by vectors of vertex coordinates.

Parameters:
x - array of x-coordinates of the triangle verteces
y - array of y-coordinates of the triangle verteces

Created by JK, 24.8.2001

Definition at line 3026 of file matrix.cpp.

double det3d ( const double *  x,
const double *  y,
const double *  z 
)

computes determinant of matrix which equals to six-multiple of tetrahedron volume which is given coordinate vectors x,y,z

The function computes determinant of the matrix

| 1 x1 y1 z1 | | 1 x2 y2 z2 | | 1 x3 y3 z3 | | 1 x4 y4 z4 |

which is equal to six-multiple of the volume of the tetrahedron whose vertices are given by coordinate vectors x,y,z.

Parameters:
x - vector of x-coordinates of tetrahedron vertices
y - vector of y-coordinates of tetrahedron vertices
z - vector of z-coordinates of tetrahedron vertices

Created by JK, 24.8.2001, Modified by JK 23.9.2011

Definition at line 3054 of file matrix.cpp.

long detm ( const matrix a,
double &  det 
)

computes determinant of the matrix a, i.e. det(A)

The function solves determinant of the matrix a by the Gauss elimination, the result is stored in the det, contents of the a is not changed

Parameters:
a is the structure of the matrix, whose determinant will being solved
det is the variable type double

Requests : a has to have the same number of the rows and columns

Return values:
0 : on succes
1 : in case incompatibility sizes of the A matrix

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2187 of file matrix.cpp.

References gause(), matrix::m, matrix::n, and print_err().

void extractblock ( double *  a,
double *  b,
long  n,
long  nr,
long  nc,
long  fri,
long  fci 
)

extracts submatrix b(nr,nc) from matrix a

The function picks up submatrix b from the original matrix a

Parameters:
a - array of original matrix of type A(n,n)
b - array of selected submatrix B(nr,nc)
n - number of columns in matrix A
nr - number of rows in the matrix B
nc - number of columns in the matrix B
fri,fci - first row and column indices of submatrix in the original matrix

Created by JK, 9.12.2004

Definition at line 3609 of file matrix.cpp.

void extractm ( matrix a,
matrix b,
long  fi,
long  ncomp 
)

extracts squared matrix b(ncomp,ncomp) from matrix a

The function extracts a squared matrix with ncomp components from matrix b and put them into matrix a.

Parameters:
a - matrix containing extracted components
b - matrix from where components are extracted
fi - first index
ncomp - number of components

Created by JK 25.2.2004

Definition at line 3584 of file matrix.cpp.

long fillcol ( double  c,
long  i,
matrix mat 
)

fills i-th column of matrix with given value c

The function fills mat's column i with the constant c.

Parameters:
c is constant, which will be used for filling entire matrix
i is index of given column
mat is the structure for the allocated matrix
Returns:
always zero

created 22.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 600 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

Referenced by mohrcoulomb::cutting_plane(), and mohrcoulomb::mc_msurf_cp().

long fillm ( double  c,
matrix mat 
)

fills contents of matrix with given value c

The function fills memory of mat's member array a with constant c.

Parameters:
c is constant, which will be used for filling entire matrix
mat is the structure for the allocated matrix
Returns:
always zero

created 22.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 552 of file matrix.cpp.

References matrix::m, and matrix::n.

long fillrow ( double  c,
long  i,
matrix mat 
)

fills i-th row of matrix with given value c

The function fills mat's row i with the constant c.

Parameters:
c is constant, which will be used for filling entire matrix
i is index of given row
mat is the structure for the allocated matrix
Returns:
always zero

created 22.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 573 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

Referenced by mohrcoulomb::cutting_plane(), and mohrcoulomb::mc_msurf_cp().

long gause ( const matrix a,
matrix b,
vector r,
vector sol,
double  zero 
)

solves system of linear algebraic equations by the Gauss elimination (A.x = r), resulting matrix is b, resulting vector is stored in r

The function solves system of linear algebraic equtions by the Gauss elimination method. The system of equation is given by matrix a, the resulting eliminated matrix is stored in b, solution is stored in the vector sol.

Parameters:
a - system matrix
b - resulting eliminated matrix
r - righthand side vector
sol - solution vector
zero - computer zero

Requests : a has to have the same number of the rows and columns b has to have same dimensions as the matrix A b has to have allocated memory array for the elements which is enough large to hold contents of the result. r has to have same dimension as number of rows in the matrix a sol has to have same dimension as the vector r sol has to have allocated memory array for the elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a and b matrices
2 : in case dimension of the matrix a < 2
3 : in case some diagonal elements is to small or zero the value which is assumed as zero is driven by macro USREPS, which is defined in the matrix.h

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1883 of file matrix.cpp.

References matrix::m, vector::n, matrix::n, and print_err().

long gause ( const matrix a,
matrix b,
vector r,
double  zero 
)

solves system of linear algebraic equations by the Gauss elimination (A.x = r), resulting matrix is b, resulting vector is stored in r

The function solves the system of equtions by the Gauss elimination. A.x = r The resulting eliminated matrix is stored b, solution vector x is stored in the vector r.

Parameters:
a - the system equation matrix
b - the resulting eliminated matrix
r - the righthand side vector/ result vector
zero - computer zero

Requests : a has to have the same number of the rows and columns b has to have same dimensions as the matrix a b has to have allocated memory array for the elements which is enough large to hold contents of the result. r has to have same dimension as number of rows in the matrix a

Return values:
0 : on succes
1 : in case incompatibility sizes of a and b matrices
2 : in case dimension of the matrix a < 2
3 : in case some diagonal elements is to small or zero the value which is assumed as zero is driven by macro USREPS, which is defined in the matrix.h

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1770 of file matrix.cpp.

References matrix::m, vector::n, matrix::n, and print_err().

long gause ( const matrix a,
matrix b,
double  zero 
)

performs Gauss elimination on the matrix a, result is stored in the matrix b

The function performs Gauss elimination on the matrix a, the result is stored in the matrix b.

Parameters:
a - elimanted matrix
b - the resulting matrix
zero - computer zero

Requests : a has to have the same number of the rows and columns b has to have same dimensions as the matrix a b has to have allocated memory array for the elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a and b matrices
2 : in case dimension of the matrix a < 2
3 : in case some diagonal elements is to small or zero the value which is assumed as zero is driven by macro USREPS, which is defined in the matrix.h

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1676 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

void gemp ( double *  a,
double *  x,
double *  y,
long  n,
long  m,
double  limit,
long  pivot 
)

solves system of algebraic equations by Gauss method with pivot selection

Function performes Gaussian elimination on matrix a. Matrix is stored as dense matrix. x and y are also dense matrices, stored by rows. right hand side vectors and vectors of unknowns are columns!

Parameters:
a - array containing matrix of the system
x - array containing vectors of unknowns
y - array containing right hand sides
n - order of matrix
m - number of right hand sides
limit - computer zero
pivot=1 - pivot is searched only when there is zero on the diagonal pivot=2 - pivot is searched every time

created 23.7.2001

Definition at line 3221 of file matrix.cpp.

References f, and print_err().

void glmatrixtransf ( matrix a,
const matrix tmat 
)

transforms tensor(matrix) a given in the global coordinate system to the local one (A_l := T^T . A_g . T)

The function transforms matrix a of element from global problem coordinate system to local element coordinate system A_l = T^T A_g T.

The transformation matrix T has to be in the form x_g = T x_l, the transormed matrix is stored in the matrix a.

Parameters:
a - matrix of element A_g (input)/resulting transformed matrix A_l(output)
tmat - transformation matrix T

Created by JK, 28.11.2006

Definition at line 2723 of file matrix.cpp.

References matrix::m, mtxm(), and mxm().

void globloc ( const double *  gv,
double *  lv,
const long *  cn,
long  n 
)

localizes components of local(element) vector to the global(problem) vector

The function localizes components of global vector gv to the local vector lv.

Parameters:
gv - array of global vector (vector of the whole problem/structure)
lv - array of local vector (vector of one element)
cn - array containing code numbers
n - number of components in the local vector

Created by JK

Definition at line 2954 of file matrix.cpp.

void glvectortransf ( const vector g,
vector l,
const matrix tmat 
)

transforms vector g given in the global coordinate system to vector l in local coordinate system (l = T^T . g)

The function transforms vector of element from global problem coordinate system to local element coordinate system l = T^T g

The transformation matrix T has to be in the form x_g = T x_l

Parameters:
g - vector in the global problem coordinate system
l - vector in the local element coordinate system
tmat - transformation matrix T

Created by JK, 28.11.2006

Definition at line 2700 of file matrix.cpp.

References matrix::m, mxv(), and tranm().

void glvectortransfblock ( vector v,
const matrix t 
)

transforms block vector v given in the global coordinate system to the local one (v_l := T^T . v_g) (number of blocks is given by dimension of T)

The function transforms vector v expressed in global coordinate system to vector expressed in local coordinate system v_l = T^T . v_g

The dimension (dim) of the problem is determined from the size of the transformation matrix T and vector v may contain several blocks with size (dim).

Parameters:
v - global vector v_g (during input), local vector v_l (during output)
t - transformation matrix T

Created by JK,

Definition at line 2821 of file matrix.cpp.

References allocv(), destrv(), vector::n, and matrix::n.

long identm ( matrix mat  ) 

sets contents of mat to identity matrix

The function sets the given matrix to the identity one..

Parameters:
mat is the structure for the allocated matrix
Returns:
always zero

Created 11.12.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 625 of file matrix.cpp.

References matrix::m, and matrix::n.

long invm ( const matrix a,
matrix b,
double  zero 
)

inverts matrix a, i.e. computes A^{-1}

The function inverts matrix A by the Gauss elimination, i.e. B := A^{-1}. The result is stored in the matrix a, contents of the a is not changed.

Parameters:
a - inverteted matrix
b - resulting matrix

Requests : a has to have the same number of the rows and columns b has to have same dimensions as the matrix a b has to have allocated memory array for the elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a and b matrices
2 : in case dimension of the matrix a < 2

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1987 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

void lgmatrixtransf ( matrix a,
const matrix tmat 
)

transforms tensor(matrix) a given in the local coordinate system to the global one (A_g := T . A_l . T^T)

The function transforms matrix a of element from local element coordinate system to global problem coordinate system A_g = T A_l T^T.

The transformation matrix T has to be in the form x_g = T x_l, the transformed matrix is stored in the matrix a

Parameters:
a - matrix of element A_l(input)/resulting transformed matrix A_g(output)
tmat - transformation matrix T

Created by JK, 28.11.2006

Definition at line 2747 of file matrix.cpp.

References matrix::m, mxm(), and mxmt().

void lgmatrixtransfblock ( matrix a,
const matrix t 
)

transforms block matrix a given in the local coordinate system to the global one (A_g := T . A_l . T^T) (number of blocks is given by dimension of T)

The function transforms matrix expressed in local coordinate system to matrix expressed in global coordinate system A_g = T . A_l . T^T

The dimension (dim) of the problem is determined from the size of the transformation matrix T. The matrix A may contain several blocks with size (dim,dim)

Parameters:
a - local matrix (during input), global matrix (during output)
t - transformation matrix

Created by JK,

Definition at line 2771 of file matrix.cpp.

References allocm(), destrm(), and matrix::n.

void lgvectortransf ( vector g,
const vector l,
const matrix tmat 
)

transforms vector l given in the local coordinate system to vector g in global coordinate system (g = T . l)

The function transforms vector of element from local element coordinate system to global problem coordinate system g = T l

The transformation matrix has to be in the form x_g = T x_l

Parameters:
g - vector in the global problem coordinate system
l - vector in the local element coordinate system
tmat - transformation matrix T

Created by JK, 28.11.2006

Definition at line 2680 of file matrix.cpp.

References mxv().

void lgvectortransfblock ( vector v,
const matrix t 
)

transforms block vector v given in the local coordinate system to the global one (v_g := T^T . v_l) (number of blocks is given by dimension of T)

The function transforms vector v expressed in local coordinate system to vector expressed in global coordinate system v_g = T * v_l

The dimension (dim) of the problem is determined from the size of the transformation matrix T and vector v may contain several blocks with size (dim)

Parameters:
v - local vector v_l (during input), global vector v_g (during output)
t - transformation matrix T

Created by JK,

Definition at line 2876 of file matrix.cpp.

References allocv(), destrv(), vector::n, and matrix::n.

void locglob ( double *  gv,
const double *  lv,
const long *  cn,
long  n 
)

localizes components of local(element) vector to the global(problem) vector

The function localizes components of local vector lv to the global vector gv.

Parameters:
gv - array of global vector (vector of the whole problem/structure)
lv - array of local vector (vector of one element)
cn - array containing code numbers
n - number of components in the local vector

Created by JK 25.7.2001

Definition at line 2928 of file matrix.cpp.

void locvecmat ( double *  mat,
const double *  vect,
const long *  cn,
long  ci,
long  m,
long  n 
)

Definition at line 2970 of file matrix.cpp.

void lu_full ( double *  a,
double *  x,
double *  y,
long  n,
double  zero,
long  tc 
)

solves system of algebraic equations by LU decomposition

function decomposes matrix A into L.U form

Parameters:
a - matrix
x - left hand side
y - right hand side
n - number of rows/columns of matrix
zero - computer zero
tc - type of computation

tc = 1 - decomposition of matrix and back substitution tc = 2 - decomposition of matrix tc = 3 - back substitution

JK, 6.1.2000

Definition at line 3398 of file matrix.cpp.

References print_err().

void mat_localize ( matrix gm,
matrix lm,
long *  rcn,
long *  ccn 
)

localizes components of local matrix lm to the global one gm

The function localizes components of local matrix lm to the global matrix gm.

Parameters:
gv - global matrix (matrix of the whole problem/structure)
lv - local matrix (vmatrix of one element)
rcn - array containing code numbers for matrix rows
ccn - array containing code numbers for matrix columns

Created by JK 25.7.2001

Definition at line 2993 of file matrix.cpp.

References matrix::m, and matrix::n.

void matassem_lsm ( double *  lsm,
vector natcoord 
)

assembles matrix of least square problem for extrapolation of nodal values from integration point values

The function assembles matrix of least square problem which is used in extrapolation of values from integration points to nodes.

Parameters:
lsm - array containing matrix
natcoord - vector containing natural coordinates of integration points

Created by JK, 10.5.2002

Definition at line 3464 of file matrix.cpp.

References vector::n.

long mixv ( const matrix a,
const vector u,
double &  vi,
long  i 
)

computes i-th component of the vector v given by product of matrix a and vector u v := A.u

The function computes i-th component of the vector v given by product of matrix a and vector u v := A.u

Parameters:
a - structure with the multiplied matrix
u - structure with the column vector
vi - i-th component of the resulting vector v
i - index of required vector component

Requests : u ,a have to be following dimensions a (m,n), u (n)

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a and vector u.

created 17.12.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1424 of file matrix.cpp.

References matrix::m, vector::n, matrix::n, and print_err().

Referenced by gl_comp_engvectortransf().

long mswapc ( matrix mat,
long  i,
long  j 
)

swaps matrix columns i and j

The function swaps mat's columns i and j

Parameters:
mat is the structure of the matrix
i number of the first column which will be swapped
j number of the second column which will be swapped
Returns:
always zero

created 22.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2360 of file matrix.cpp.

References matrix::m.

long mswapr ( matrix mat,
long  i,
long  j 
)

swaps matrix rows i and j

The function swaps mat's rows i and j

Parameters:
mat is the structure of the matrix
i number of the first row which will be swapped
j number of the second row which will be swapped
Returns:
always zero

created 22.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2386 of file matrix.cpp.

References matrix::n.

void mtxm ( const double *  a,
const double *  b,
double *  c,
long  l,
long  m,
long  n 
)

multiplies two matrices C := A^T.B given by arrays

The function multiplies transposed matrix given by 1D double array a from right with matrix given by 1D double array b, the result is stored in c A[l,m]^T . B[m,n] = C[l,n]

Parameters:
a is the 1D double array of first matrix, from which is multiplied
b is the 1D double array of the second=multiplicating matrix
c is the 1D double array of the result matrix

Requests : a, b and c have to be following dimensions a (l*m), b (m*n), c (l*n)

Return values:
0 : always

created 13.4.1998 by JK

Definition at line 1044 of file matrix.cpp.

long mtxm ( const matrix a,
const matrix b,
matrix c 
)

multiplies two matrices C := A^T.B

The function multiplies transposed matrix given by a from right with matrix given by b, the result is stored in c

Parameters:
a is the structure of the matrix, from which is multiplied
b is the structure of the multiplicating matrix
c is the structure of the result matrix

Requests : a, b and c have to be following dimensions a (n,m), b (n,p), c (m,p) c has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a, b and c matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1004 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

void mtxmccr ( const double *  a,
const double *  b,
double *  c,
long  l,
long  m,
long  n 
)

multiplies transposed matrix a from right by matrix b, C := A^T . B, matrices are stored by columns

The function computes matrix product C := A^T . B, all matrices are stored by columns.

Parameters:
a - array with matrix A(m,n)
b - array with matrix B(m,p)
c - array with matrix C(n,p)
l -
m -
n -

Created by JK, 7.12.1998

Definition at line 1077 of file matrix.cpp.

void mtxv ( const double *  a,
const double *  b,
double *  c,
long  m,
long  n 
)

multiplies transposed matrix a given by array from right by vector u given by array v := A^T.u

The function multiplies transposed matrix by vector from right, i.e. c:= A^T.b.

Parameters:
a(m,n) - matrix A given by array
b(m,1) - vector b given by array
c(n,1) - resulting vector c given by array
m - number of matrix rows
n - number of matrix columns, number of vector components

Created by JK, 19.2.1997

Definition at line 1495 of file matrix.cpp.

long mtxv ( const matrix a,
const vector u,
vector v 
)

multiplies transposed matrix a from right by vector u v := A^T.u

The function multiplies transposed matrix given by a with column vector u, the result is stored in column vector v.

Parameters:
a is the structure of the matrix, which is multiplied
u is the structure of the column vector
v is the structure of the result column vector

Requests : u ,a and v have to be following dimensions u (m), a (n,m), c (n) v has to have allocated memory array for the elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1463 of file matrix.cpp.

References matrix::m, matrix::n, vector::n, and print_err().

void mtxvc ( const double *  a,
const double *  b,
double *  c,
long  m,
long  n 
)

The function multiplies transposed matrix a by the vector from right, i.e. c:= A^T.b, The matrix a is stored by columns.

Parameters:
a - array with matrix a(m,n) stored by columns
b - array with vector b(m)
m - number of matrix rows / vector b components
n - number of matrix columns / vector c components
c - array with vector c(n)

Created 7.12.1998 by JK

Definition at line 1524 of file matrix.cpp.

void mxm ( const double *  a,
const double *  b,
double *  c,
long  l,
long  m,
long  n 
)

multiplies two matrices C := A.B given by arrays

The function multiplies matrix given by 1D double array a from right with matrix given by 1D double array b, the result is stored in 1D double array c. A[l,m] . B[m,n] = C[l,n]

Parameters:
a is the 1D double array of first matrix, from which is multiplied
b is the 1D double array of the second=multiplicating matrix
c is the 1D double array of the result matrix

Requests : a, b and c have to be following dimensions a (l*m), b (m*n), c (l*n)

Return values:
0 : always

created 19.2.1997 by JK

Definition at line 880 of file matrix.cpp.

long mxm ( const matrix a,
const matrix b,
matrix c 
)

multiplies two matrices C := A.B

The function multiplies matrix given by a from right with matrix given by b, the result is stored in c

Parameters:
a is the structure of the matrix, from which is multiplied
b is the structure of the multiplicating matrix
c is the structure of the result matrix

Requests : a, b and c have to be following dimensions a (m,n), b (n,p), c (m,p) c has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a, b and c matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 840 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

void mxmt ( const double *  a,
const double *  b,
double *  c,
long  l,
long  m,
long  n 
)

multiplies two matrices C := A.B^T given by arrays

The function multiplies matrix given by 1D double array a from right with transposed matrix given by 1D double array b, the result is stored in c A[l,m] . B[m,n]^T = C[l,n]

Parameters:
a is the 1D double array of first matrix, from which is multiplied
b is the 1D double array of the second=multiplicating matrix
c is the 1D double array of the result matrix

Requests : a, b and c have to be following dimensions a (l*m), b (m*n), c (l*n)

Return values:
0 : always

created 9.11.1999 by JK

Definition at line 964 of file matrix.cpp.

long mxmt ( const matrix a,
const matrix b,
matrix c 
)

multiplies two matrices C := A.B^T

The function multiplies matrix given by a from right with transposed matrix given by b, the result is stored in c

Parameters:
a is the structure of the matrix, from which is multiplied
b is the structure of the multiplicating matrix
c is the structure of the result matrix

Requests : a, b and c have to be following dimensions a (m,n), b (p,n), c (m,p) c has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a, b and c matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 924 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

void mxv ( const double *  a,
const double *  b,
double *  c,
long  m,
long  n 
)

multiplies matrix a given by array from right by vector u given by array v := A.u

funkce nasobi matici A(m,n) s vektorem b(n,1) A(m,n).b(n,1)=c(m,1)

Parameters:
c - vysledny vektor
a - matice A
b - vektor b
m 
n - rozmer matice A

created 19.2.1997

Definition at line 1359 of file matrix.cpp.

long mxv ( const matrix a,
const vector u,
vector v 
)

multiplies matrix a from right by vector u v := A.u

The function multiplies matrix given by a with column vector u, the result is stored in column vector v.

Parameters:
a is the structure of the matrix, which is multiplied
u is the structure of the column vector
v is the structure of the result column vector

Requests : u ,a and v have to be following dimensions v (m), a (m,n), u (n) v has to have allocated memory array for the elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1326 of file matrix.cpp.

References matrix::m, vector::n, matrix::n, and print_err().

void mxvc ( const double *  a,
const double *  b,
double *  c,
long  m,
long  n 
)

multiplies matrix a stored by columns from the right by vector b c := A.b

The function multiplies matrix a stored by columns in the array by the vector b stored in array. The resulting vector is stored in the array c. (c := A.b)

Parameters:
a - array with matrix components stored by columns
b - array with components of vector b
m - the number of matrix rows
n - the number of matrix/vector columns
c - array of the resulting vector

Created 15.12.1998 by Jaroslav Kruis

Definition at line 1390 of file matrix.cpp.

void nnj ( double *  m,
const double *  n,
double  jac,
long  mm,
long  nn 
)

computes matrix product M += N^T . N . jac

The function evaluates matrix product N^T . N . jac, the result is added to the matrix m.

Parameters:
m - array of resulting matrix M += N^T . N. jac
n - array of matrix N
jac - jacobian (constant)
mm - number of rows of matrix n
nn - number of columns of matrix n

Created by JK 9.7.2001

Definition at line 2614 of file matrix.cpp.

References ll.

void nnjac ( matrix d,
const matrix a,
const matrix b,
double  jac 
)

computes matrix product D += A^T . B . jac

The function evaluates matrix product A^T . B . jac, the result is added to the matrix d.

Parameters:
d - resulting matrix D += A^T . B . jac
a - matrix A
b - matrix B
jac - jacobian (constant)

Created by JK 10.5.2002

Definition at line 2655 of file matrix.cpp.

References addm(), allocm(), cmulm(), destrm(), mtxm(), and matrix::n.

void princ_val ( matrix v,
vector pval,
matrix pvect,
long  ni,
double  err,
double  zero,
long  n,
long  normalize 
)

computes eigenvalues(principal values) of the given matrix v by Jacobi method

Function computes principal values and vectors of tensors of 2D and 3D problems .

Parameters:
v - matrix of second or third order tensor
pval - vector of sorted principal values (from min to max)
pvect - matrix of principal direction vectors stored in columns
ni - number of iterations in jacobi method
err - required error of Jacobi method
zero - computer zero for explicit computation method
n - order of problem (only 2 for 2D problem and 3 for 3D problem are supported)
normalize - switch for normalizing of matrix in Jacobi method
Returns:
eigenvectors create columns of the matrix pvect and eigenvalues are stored and sorted in pval vector

Created 29.8.2001 Modified 17.8.2007 TKo (added normalizing in Jacobi method)

Definition at line 3108 of file matrix.cpp.

References vector::a, jacobi(), mswapc(), and print_err().

long printm ( FILE *  out,
const gfmatrix a 
)

prints the gfmatrix to the file

The function prints out the contents of the matrix a to the file given by out.

Parameters:
out is the structure with opened file for matrix output
a is the structure of the matrix which is printed
Returns:
always zero.

Created 10.12.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2332 of file matrix.cpp.

References gfmatrix::m, and gfmatrix::n.

long printm ( FILE *  out,
const matrix a 
)

prints the matrix to the file

The function prints out the contents of the matrix a to the file given by out.

Parameters:
out is the structure with opened file for matrix output
a is the structure of the matrix which is printed
Returns:
always zero.

created 10.12.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2307 of file matrix.cpp.

References matrix::m, and matrix::n.

long printm ( const matrix a,
FILE *  out,
int  prec,
int  width 
)

prints the matrix to the file with the given width and precision

The function prints out the contents of the matrix a to the stdout file, in precision 3 digits, width of the number is 11 optionally, can be specified file to print, precision and the field width.

Parameters:
a is the structure of the matrix which is printed

Optionally :

Parameters:
out is the structure with opened file for matrix output (default is stdout)
prec is desired precision of the matrix elements (default is 3)
width is desired width of the matrix elements on the output (default is 11)
Returns:
always zero.

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2282 of file matrix.cpp.

References matrix::m, and matrix::n.

long readm ( XFILE in,
gfmatrix a 
)

reads gfmatrix from the file in

The function reads the contents of the gfmatrix a from the opened text file.

Parameters:
in - pointer to the opened XFILE
a - is the structure of the matrix which is read
Returns:
always zero.

Created 10.12.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2251 of file matrix.cpp.

References gfmatrix::m, and gfmatrix::n.

long readm ( XFILE in,
matrix a 
)

reads matrix from the file in

The function reads the contents of the matrix a from the opened text file.

Parameters:
in - pointer to the opened XFILE
a - is the structure of the matrix which is read
Returns:
always zero.

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 13.10.2014, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 2227 of file matrix.cpp.

References matrix::m, matrix::n, and xfscanf().

long reallocm ( long  m,
long  n,
imatrix mat 
)

reallocates imatrix to the dimensions m x n

The function reallocates memory for the member array a. If the memory requirements are less or equal to the allocated one, the memory is not reallocated and the matrix dimensions are changed only.

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated imatrix

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 486 of file matrix.cpp.

References imatrix::a, imatrix::m, memset(), imatrix::n, and print_err().

long reallocm ( long  m,
long  n,
matrix mat 
)

reallocates matrix to the dimensions m x n

The function reallocates memory for the member array a. If the memory requirements are less or equal to the allocated one, the memory is not reallocated and the matrix dimensions are changed only.

Parameters:
m is the number of rows
n is the number of columns
mat is the structure for allocated matrix

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 448 of file matrix.cpp.

References matrix::a, matrix::m, memset(), matrix::n, and print_err().

void rhsassem_lsm ( double *  rhs,
vector natcoord,
vector values 
)

assembles right hand side vector of least square problem for extrapolation of nodal values from integration point values

The function assembles right hand side of least square problems. It is used in extrapolation of values from integration points to nodes.

Parameters:
rhs - array containing matrix
natcoord - vector of natural coordinates of integration points
values - vector of values in integration point

Created by JK, 10.5.2002

Definition at line 3495 of file matrix.cpp.

References vector::n.

void solve_lsm ( double *  lsm,
double *  lhs,
double *  rhs,
double  zero,
long  n,
long  m 
)

solves system of equation given by least square problem for extrapolation of nodal values from integration point values

The function solves least square problems. It is used in extrapolation of values from integration points to nodes

Parameters:
lsm - array of matrix of least square problem
lhs - array of left hand side (contains unknown coefficients)
rhs - array of right hand side
zero - computer zero
n - number of rows/columns of matrix
m - number of values

Created by JK, 10.5.2002

Definition at line 3536 of file matrix.cpp.

References lu_full().

long subm ( const matrix a,
const matrix b,
matrix c 
)

substracts 2 matrices C := A-B

The function subtracts matrix given by b from matrix given by a, the result is stored in the c

Parameters:
a is the structure of the matrix, from which is subtracted
b is the structure of subtracted matrix
c is the structure of the result matrix

Requests : a, b and c have to be same dimensions, c has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a, b and c matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 765 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

long tensprd ( const vector a,
const vector b,
matrix c 
)

The function performs tensor product of two vectors.

Parameters:
a is the structure of the first vector
b is the structure of the second vector
c is the structure of the result matrix

Requests : c has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a, b and c matrices

created

Definition at line 801 of file matrix.cpp.

References matrix::m, matrix::n, vector::n, and print_err().

long tranm ( matrix a  ) 

transposes matrix A := A^T

The function transposes matrix given by a, the result is stored in a.

Parameters:
a is the structure of the matrix, from which is transposed
Return values:
always zero

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1205 of file matrix.cpp.

References matrix::a, matrix::m, and matrix::n.

long tranm ( const matrix a,
matrix at 
)

transposes matrix A_t := A^T

The function transposes matrix given by a, the result is stored in at

Parameters:
a is the structure of the matrix, from which is transposed
at is the structure of the result matrix

Requests : a and at have to be same dimensions, at has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of a and at matrices

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1175 of file matrix.cpp.

References matrix::m, matrix::n, and print_err().

long vxm ( const vector u,
const matrix a,
vector v 
)

multiplies matrix a from left by row vector u v := u.A

The function multiplies row vector given by u with matrix given by a, the result is stored in row vector v

Parameters:
u is the structure of the row vector
a is the structure of the matrix, which is multiplied
v is the structure of the result row vector

Requests : u ,a and v have to be following dimensions u (m), a (m,n), c (n) v has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1243 of file matrix.cpp.

References matrix::m, matrix::n, vector::n, and print_err().

long vxmt ( const vector u,
const matrix a,
vector v 
)

multiplies transposed matrix a from left by row vector u v := u.A^T

The function multiplies row vector given by u with transposed matrix given by a, the result is stored in row vector v.

Parameters:
u is the structure of the row vector
a is the structure of the matrix, which is multiplied
v is the structure of the result row vector

Requests : u ,a and v have to be following dimensions u (m), a (n,m), c (n) v has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1284 of file matrix.cpp.

References matrix::m, matrix::n, vector::n, and print_err().

long vxmxv ( const double *  v,
const double *  m,
long  dim,
double &  answer 
)

multiplies row vector v given by array with symetric matrix m given by array and column vector v, i.e. answer := (v^T).A.v

The function multiplies row vector given by "v" with symetric matrix given by "m" with column vector given by "v" -> v.m.v , the result is stored in answer.

Parameters:
v is pointer on double array
m is pointer on double array
dim is size of "v" & "m"
answer is the real number

Requests : "v" and "m" have to be following dimensions: v (dim), m (dim,dim)

Return values:
0 : always

created 25.10.2002, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 1635 of file matrix.cpp.

long vxmxv ( const vector v,
const matrix m,
double &  answer 
)

multiplies row vector v with symetric matrix m and column vector v, i.e. answer := (v^T).A.v

The function multiplies row vector given by "v" with symetric matrix given by "m" with column vector given by "v" -> v.m.v , the result is stored in answer.

Parameters:
v is the structure of the row vector
m is the structure of the matrix
answer is the real number

Requests : "v" and "m" have to be following dimensions: v (n), m (n,n)

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a

created 25.10.2002, Ladislav Svoboda, termit@cml.fsv.cvut.cz

Definition at line 1595 of file matrix.cpp.

References matrix::m, matrix::n, vector::n, and print_err().

long vxv ( const vector u,
const vector v,
matrix a 
)

multiplies column vector u with row vector, i.e. it performs tensor product, A := u * v

The function multiplies column vector given by u with row vector given by v, the result is stored in matrix a.

Parameters:
u is the structure of the column vector
v is the structure of the row vector
a is the structure of the result matrix

Requests : u ,a and v have to be following dimensions u (m), v(n), a(m,n) a has to have allocated memory array for elements which is enough large to hold contents of the result.

Return values:
0 : on succes
1 : in case incompatibility sizes of the matrix a

created 29.5.2000, Tomas Koudelka, koudelka@cml.fsv.cvut.cz modified 9.5.2001, Tomas Koudelka, koudelka@cml.fsv.cvut.cz

Definition at line 1560 of file matrix.cpp.

References matrix::m, matrix::n, vector::n, and print_err().


Generated by  doxygen 1.6.2