Elmer FEM solver
Elmer is an open source finite element software for multiphysical problems
 All Classes Files Functions Variables Typedefs Macros Groups Pages
crsmatrix Module Reference

Public Member Functions

integer function crs_search (N, Array, VALUE)
 
type(matrix_t) function, pointer crs_transpose (A)
 
subroutine crs_mergematrix (A, B, PermA, PermB)
 
subroutine crs_zeromatrix (A)
 
subroutine crs_zerorow (A, n)
 
subroutine crs_sortmatrix (A, ValuesToo)
 
subroutine crs_sortbasicmatrix (A, ValuesToo)
 
subroutine crs_makematrixindex (A, i, j)
 
subroutine crs_addtomatrixelement (A, i, j, VALUE)
 
subroutine crs_setmatrixelement (A, i, j, VALUE)
 
real(kind=dp) function crs_getmatrixelement (A, i, j)
 
real(kind=dp) function crs_changematrixelement (A, i, j, NewValue)
 
subroutine crs_moverow (A, n1, n2, coeff)
 
subroutine crs_gluelocalmatrix (A, N, Dofs, Indeces, LocalMatrix)
 
subroutine crs_gluelocalsubmatrix (A, row0, col0, Nrow, Ncol, RowInds, ColInds, RowDofs, ColDofs, LocalMatrix)
 
subroutine crs_setsymmdirichlet (A, b, n, val)
 
real(kind=dp) function crs_rowsum (A, k)
 
subroutine crs_rowsuminfo (A, Values)
 
type(matrix_t) function, pointer crs_creatematrix (N, Total, RowNonzeros, Ndeg, Reorder, AllocValues)
 
subroutine crs_printmatrix (A)
 
subroutine crs_printrhs (A)
 
subroutine crs_matrixvectormultiply (A, u, v)
 
subroutine crs_transposematrixvectormultiply (A, u, v)
 
subroutine crs_complexmatrixvectormultiply (A, u, v)
 
subroutine crs_applyprojector (PMatrix, u, uperm, v, vperm, Trans)
 
subroutine crs_diagprecondition (u, v, ipar)
 
subroutine crs_complexdiagprecondition (u, v, ipar)
 
subroutine crs_blockdiagonal (A, B, Blocks)
 
subroutine crs_removezeros (A)
 
subroutine crs_fctloworder (A)
 
subroutine crs_copymatrixtopology (A, B)
 
subroutine crs_blockmatrixpick (A, B, Blocks, Nrow, Ncol)
 
subroutine crs_blockmatrixpick2 (A, B, BlockStruct, Nrow, Ncol)
 
logical function crs_copymatrixprec (A, B)
 
logical function crs_incompletelu (A, ILUn)
 
logical function crs_complexincompletelu (A, ILUn)
 
logical function crs_ilut (A, TOL)
 
logical function crs_complexilut (A, TOL)
 
subroutine crs_luprecondition (u, v, ipar)
 
subroutine crs_complexluprecondition (u, v, ipar)
 
subroutine crs_lusolve (N, A, b)
 
subroutine crs_complexlusolve (N, A, b)
 
subroutine crs_matrixvectorprod (u, v, ipar)
 
subroutine crs_complexmatrixvectorprod (u, v, ipar)
 
subroutine crs_inspectmatrix (A)
 

Member Function/Subroutine Documentation

subroutine crsmatrix::crs_addtomatrixelement ( type(matrix_t)  A,
integer, intent(in)  i,
integer, intent(in)  j,
real(kind=dp), intent(in)  VALUE 
)

Add a given value to an element of a CRS format matrix.

Parameters
aStructure holding the matrix
[in]irow number of the matrix element
[in]jcolumn number of the matrix element
[in]valuevalue to be added to the matrix element

References crs_search().

Referenced by solverutils::addtomatrixelement(), crs_moverow(), gebhardtfactors(), nodaldisplacementpenalty(), radiationfactors(), and vankacreate().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_applyprojector ( type(matrix_t)  PMatrix,
real(kind=dp), dimension(:)  u,
integer, dimension(:), pointer  uperm,
real(kind=dp), dimension(:)  v,
integer, dimension(:), pointer  vperm,
logical, optional  Trans 
)

Referenced by applyprojector(), and gmgsweep().

Here is the caller graph for this function:

subroutine crsmatrix::crs_blockdiagonal ( type(matrix_t), intent(in)  A,
type(matrix_t)  B,
integer, intent(in)  Blocks 
)

Pics the block diagonal entries from matrix A to build matrix B.

Parameters
[in]aThe initial matrix
bThe block diagonal matrix
[in]blocksNumber of blocks used in the decomposition

Referenced by itersolve::itersolver().

Here is the caller graph for this function:

subroutine crsmatrix::crs_blockmatrixpick ( type(matrix_t), intent(in)  A,
type(matrix_t)  B,
integer, intent(in)  Blocks,
integer, intent(in)  Nrow,
integer, intent(in)  Ncol 
)

Pics a block from matrix A to build matrix B. It is assumed that the matrix is split into given number of equally sized blocks.

Parameters
[in]aInitial matrix
bSubmatrix picked from the larger matrix
[in]blocksNumber of blocks in the initial matrix
[in]nrowRow to be picked
[in]ncolColumn to be picked

References messages::fatal(), and messages::warn().

Referenced by blocksolve::blockpickmatrix().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_blockmatrixpick2 ( type(matrix_t), intent(in)  A,
type(matrix_t)  B,
integer, dimension(:), pointer  BlockStruct,
integer, intent(in)  Nrow,
integer, intent(in)  Ncol 
)

Pics a block from matrix A to build matrix B. This subroutine enables the use of nontrivial block decompositions.

Parameters
[in]aInitial matrix
bSubmatrix picked from the larger matrix
blockstructBlock decomposition structure of the initial matrix
[in]nrowRow to be picked
[in]ncolColumn to be picked

References messages::fatal(), and messages::warn().

Referenced by blocksolve::blockpickmatrix().

Here is the call graph for this function:

Here is the caller graph for this function:

real(kind=dp) function crsmatrix::crs_changematrixelement ( type(matrix_t), intent(in)  A,
integer, intent(in)  i,
integer, intent(in)  j,
real(kind=dp), intent(in)  NewValue 
)

Get a given matrix entry from CRS format matrix and replace it with a new value.

Parameters
[in]aStructure holding the matrix
[in]irow number of the matrix element
[in]jcolumn number of the matrix element
[in]newvalueValue to be set

References crs_search().

Referenced by solverutils::changematrixelement().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_complexdiagprecondition ( complex(kind=dp), dimension(*), intent(out)  u,
complex(kind=dp), dimension(*), intent(in)  v,
integer, dimension(*)  ipar 
)

Diagonal preconditioning of a CRS format matrix for complex valued matrix equations. Matrix is accessed from a global variable GlobalMatrix.

Parameters
[out]uResulting approximate solution after preconditioning
[in]vGiven right-hand-side
iparstructure holding info from (HUTIter-iterative solver package)

References generalutils::sortf().

Referenced by itersolve::itersolver().

Here is the call graph for this function:

Here is the caller graph for this function:

logical function crsmatrix::crs_complexilut ( type(matrix_t)  A,
real(kind=dp), intent(in)  TOL 
)

Buids an incomplete (ILUT) factorization for an iterative solver preconditioner. Complex matrix version.

Parameters
aStructure holding input matrix, will also hold the factorization on exit.
[in]tolDrop toleranece: if ILUT(i,j) <= NORM(A(i,:))*TOL the value is dropped.

References complexcomputeilut(), and messages::info().

Referenced by eigensolve::arpackeigensolvecomplex(), and itersolve::itersolver().

Here is the call graph for this function:

Here is the caller graph for this function:

logical function crsmatrix::crs_complexincompletelu ( type(matrix_t)  A,
integer, intent(in)  ILUn 
)

Buids an incomplete (ILU(n)) factorization for an iterative solver preconditioner. Complex matrix version.

Parameters
aStrcture holding input matrix, will also hold the factorization on exit.
[in]ilunOrder of fills allowed 0-9

References messages::fatal(), messages::info(), and initializecomplexilu1().

Referenced by bicgstablouteriteration(), bicgstabouteriteration(), gcrouteriteration(), innerouteriteration(), itersolve::itersolver(), and monolithicsolve().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_complexluprecondition ( complex(kind=dp), dimension(huti_ndim), intent(out)  u,
complex(kind=dp), dimension(huti_ndim), intent(in)  v,
integer, dimension(*), intent(in)  ipar 
)

Incomplete factorization preconditioner solver for a CRS format matrix. Matrix is accessed from a global variable GlobalMatrix. The Complex version.

Parameters
[in]iparstructure holding info from (HUTIter-iterative solver package)
[in]vRight-hand-side vector
[out]uSolution vector

References crs_complexlusolve().

Referenced by itersolve::itersolver().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_complexlusolve ( integer, intent(in)  N,
type(matrix_t), intent(in)  A,
complex(kind=dp), dimension(n)  b 
)

Solve a complex system Ax=b after factorization A=LUD has been done. This routine is meant as a part of a preconditioner for an iterative solver.

Parameters
[in]aStructure holding input matrix
[in]nSize of the system
bon entry the RHS vector, on exit the solution vector.

References complexlusolve().

Referenced by eigensolve::arpackeigensolvecomplex(), ccg(), and crs_complexluprecondition().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_complexmatrixvectormultiply ( type(matrix_t), intent(in)  A,
complex(kind=dp), dimension(*), intent(in)  u,
complex(kind=dp), dimension(*), intent(out)  v 
)

Matrix vector product (v = Au) for a matrix given in CRS format assuming complex valued matrix equation.

Parameters
[in]uVector to be multiplied
[out]vResult vector
[in]aStructure holding matrix

Referenced by amgsweep(), eigensolve::arpackeigensolvecomplex(), cmgsweep(), and mgcmv().

Here is the caller graph for this function:

subroutine crsmatrix::crs_complexmatrixvectorprod ( complex(kind=dp), dimension(huti_ndim), intent(in)  u,
complex(kind=dp), dimension(huti_ndim)  v,
integer, dimension(*), intent(in)  ipar 
)

Complex matrix vector product (v = Au) for a matrix given in CRS format. The matrix is accessed from a global variable GlobalMatrix.

Parameters
[in]iparStructure holding info HUTIter-iterative solver package
[in]uvector to multiply u
vresult vector

References complexmatvec().

Referenced by itersolve::itersolver().

Here is the call graph for this function:

Here is the caller graph for this function:

logical function crsmatrix::crs_copymatrixprec ( type(matrix_t), intent(in)  A,
type(matrix_t)  B 
)

Copies some preconditioning structures from matrix A to B. The intent is to allow saving of memory and CPU time for cases where similar preconditioning could be used for many components.

Parameters
[in]aInitial matrix
bNew matrix with the same preconditioner

References messages::info().

Referenced by blocksolve::blockmatrixprec(), and blocksolve::blockstandarditer().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_copymatrixtopology ( type(matrix_t), intent(in)  A,
type(matrix_t)  B 
)

Copies the matrix topology from matrix A to build matrix B. Note that the topology is really reused, so that if matrix A is destroyed also matrix B becomes unusable.

Parameters
[in]aInitial matrix
bNew matrix with same matrix topology

References messages::fatal(), and messages::info().

Referenced by blocksolve::blockpickmatrix(), and blocksolve::blockprecmatrix().

Here is the call graph for this function:

Here is the caller graph for this function:

type(matrix_t) function, pointer crsmatrix::crs_creatematrix ( integer, intent(in)  N,
integer, intent(in)  Total,
integer, dimension(:), intent(in)  RowNonzeros,
integer, intent(in)  Ndeg,
integer, dimension(:), intent(in)  Reorder,
logical, intent(in)  AllocValues 
)

Create the structures required for a CRS format matrix.

Parameters
[in]nNumber of rows in the matrix
[in]totalTotal number of nonzero entries in the matrix
[in]ndegNegrees of freedom
[in]rownonzerosNumber of nonzero entries in rows of the matrix
[in]reorderPermutation index for bandwidth reduction
[in]allocvaluesShould the values arrays be allocated ?

References generalutils::allocatematrix(), and messages::fatal().

Referenced by elementutils::creatematrix(), gebhardtfactors(), and radiationfactors().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_diagprecondition ( real(kind=dp), dimension(*), intent(out)  u,
real(kind=dp), dimension(*), intent(in)  v,
integer, dimension(*)  ipar 
)

Diagonal preconditioning of a CRS format matrix. Matrix is accessed from a global variable GlobalMatrix. Note that if the matrix has been scaled so that the diagonal entries are already ones, this subroutine is obsolite.

Parameters
[out]uResulting approximate solution after preconditioning
[in]vGiven right-hand-side
iparstructure holding info from (HUTIter-iterative solver package)

References generalutils::sortf().

Referenced by itersolve::itersolver().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_fctloworder ( type(matrix_t)  A)

Makes a algebraic lower order scheme assuming steady state advection-diffusion equation. This can be applied together with flux corrected transport (FCT) scheme. Also creates a lumped mass to MassValuesLumped and saves the original stiffness matrix values to BulkValues.

For more information see, for example, Dmitri Kuzmin (2008): "Explicit and implicit FEM-FCT algorithms with flux linearization"

Parameters
aInitial higher order matrix

References crs_rowsuminfo(), messages::fatal(), and messages::info().

Referenced by defutils::defaultfinishassembly().

Here is the call graph for this function:

Here is the caller graph for this function:

real(kind=dp) function crsmatrix::crs_getmatrixelement ( type(matrix_t), intent(in)  A,
integer, intent(in)  i,
integer, intent(in)  j 
)

Get a given matrix entry from CRS format matrix.

Parameters
[in]aStructure holding the matrix
[in]irow number of the matrix element
[in]jcolumn number of the matrix element

References crs_search().

Referenced by solverutils::getmatrixelement(), vanka(), and vankacreate().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_gluelocalmatrix ( type(matrix_t)  A,
integer, intent(in)  N,
integer, intent(in)  Dofs,
integer, dimension(:), intent(in)  Indeces,
real(kind=dp), dimension(:,:), intent(in)  LocalMatrix 
)

Add a set of values (.i.e. element stiffness matrix) to a CRS format matrix.

Parameters
aStructure holding matrix
[in]localmatrixA (N x Dofs) x ( N x Dofs) matrix holding the values to be added to the CRS format matrix
[in]nNumber of nodes in element
[in]dofsNumber of degrees of freedom for one node
[in]indecesMaps element node numbers to global (or partition) node numbers

References localmatrix().

Referenced by modeldescription::getnodalelementsize(), freesurface::poissonsolve(), solverutils::updateglobalequations(), updateglobalpreconditioner(), and solverutils::updatemassmatrix().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_gluelocalsubmatrix ( type(matrix_t)  A,
integer, intent(in)  row0,
integer, intent(in)  col0,
integer, intent(in)  Nrow,
integer, intent(in)  Ncol,
integer, dimension(:), intent(in)  RowInds,
integer, dimension(:), intent(in)  ColInds,
integer, intent(in)  RowDofs,
integer, intent(in)  ColDofs,
real(kind=dp), dimension(:,:), intent(in)  LocalMatrix 
)

Add a set of values (.i.e. element stiffness matrix) to a CRS format matrix. For this matrix the entries are ordered so that first for one dof you got all nodes, and then for second etc. There may be on offset to the entries making the subroutine suitable for coupled monolithic matrix assembly.

Parameters
[in]localmatrixA (Nrow x RowDofs) x ( Ncol x ColDofs) matrix holding the values to be added to the CRS format matrix
aStructure holding matrix
[in]nrowNumber of nodes for the row variable
[in]ncolNumber of nodes for the column variable
[in]rowdofsNumber of dofs for row variable
[in]coldofsNumber of dofs for column variable
[in]col0Offset for column variable
[in]row0Offset for row variable
[in]rowindsPermutation of the row dofs
[in]colindsPermutation of the column dofs

References localmatrix().

Referenced by solverutils::gluelocalsubmatrix().

Here is the call graph for this function:

Here is the caller graph for this function:

logical function crsmatrix::crs_ilut ( type(matrix_t)  A,
real(kind=dp), intent(in)  TOL 
)

Buids an incomplete (ILUT) factorization for an iterative solver preconditioner. Real matrix version.

Parameters
aStructure holding input matrix, will also hold the factorization on exit.
[in]tolDrop toleranece: if ILUT(i,j) <= NORM(A(i,:))*TOL the value is dropped.

References computeilut(), and messages::info().

Referenced by multigrid::amgsolve(), eigensolve::arpackeigensolve(), eigensolve::arpackstabeigensolve(), multigrid::cmgsolve(), multigrid::gmgsolve(), itersolve::itersolver(), paralleleigensolve::parallelarpackeigensolve(), and multigrid::pmgsolve().

Here is the call graph for this function:

Here is the caller graph for this function:

logical function crsmatrix::crs_incompletelu ( type(matrix_t)  A,
integer, intent(in)  ILUn 
)

Builds an incomplete (ILU(n)) factorization for a iterative solver preconditioner. Real matrix version.

Parameters
aStrcture holding input matrix, will also hold the factorization on exit.
[in]ilunOrder of fills allowed 0-9

References messages::fatal(), messages::info(), initializeilu1(), and messages::warn().

Referenced by multigrid::amgsolve(), eigensolve::arpackdampedeigensolve(), multigrid::cmgsolve(), multigrid::gmgsolve(), itersolve::itersolver(), paralleleigensolve::parallelarpackeigensolve(), and multigrid::pmgsolve().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_inspectmatrix ( type(matrix_t)  A)

Check the matrix for correctness.

Parameters
aStructure holding the matrix
subroutine crsmatrix::crs_luprecondition ( real(kind=dp), dimension(huti_ndim), intent(out)  u,
real(kind=dp), dimension(huti_ndim), intent(in)  v,
integer, dimension(*), intent(in)  ipar 
)

Incomplete factorization preconditioner solver for a CRS format matrix. Matrix is accessed from a global variable GlobalMatrix. The real version.

Parameters
[in]iparstructure holding info from (HUTIter-iterative solver package)
[in]vRight-hand-side vector
[out]uSolution vector

References crs_lusolve().

Referenced by itersolve::itersolver().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_lusolve ( integer, intent(in)  N,
type(matrix_t), intent(in)  A,
double precision, dimension(n)  b 
)

Solve a system (Ax=b) after factorization A=LUD has been done. This routine is meant as a part of a preconditioner for an iterative solver.

Parameters
[in]aStructure holding input matrix
[in]nSize of the system
bon entry the RHS vector, on exit the solution vector.

References lusolve().

Referenced by bicg(), paralleleigensolve::bicgpareigen(), bicgstab(), cg(), paralleleigensolve::cgpareigen(), crs_luprecondition(), and eigensolve::eigenbicg().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_makematrixindex ( type(matrix_t)  A,
integer, intent(in)  i,
integer, intent(in)  j 
)

Fill in the column number to a CRS format matrix (values are not affected in any way).

Parameters
aStructure holding matrix
[in]irow number of the matrix element
[in]jcolumn number of the matrix element

References messages::error().

Referenced by elementutils::creatematrix(), crs_matrixmatrixmultiply(), gebhardtfactors(), and radiationfactors().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_matrixvectormultiply ( type(matrix_t), intent(in)  A,
real(kind=dp), dimension(*), intent(in)  u,
real(kind=dp), dimension(*), intent(out)  v 
)

Matrix vector product (v = Au) for a matrix given in CRS format.

Parameters
[in]uVector to be multiplied
[out]vResult vector
[in]aStructure holding matrix

Referenced by amgsweep(), eigensolve::arpackeigensolve(), eigensolve::arpackstabeigensolve(), blocksolve::blockmatrixprec(), blocksolve::blockmatrixvectorprod(), blocksolve::blockupdaterhs(), eigensolve::checkresiduals(), cmgsweep(), eigensolve::eigenmgmv1(), eigensolve::eigenmgmv2(), eliminateperiodic(), freesurfacesolver(), solverutils::matrixvectormultiply(), mgmv(), mymv(), pmgsweep(), pmv(), preconditioningiteration(), rigidbody(), stokessolver(), and itersolve::stopc().

Here is the caller graph for this function:

subroutine crsmatrix::crs_matrixvectorprod ( real(kind=dp), dimension(huti_ndim), intent(in)  u,
real(kind=dp), dimension(huti_ndim)  v,
integer, dimension(*), intent(in)  ipar 
)

Matrix vector product (v = Au) for a matrix given in CRS format. The matrix is accessed from a global variable GlobalMatrix.

Parameters
[in]iparStructure holding info HUTIter-iterative solver package
[in]uvector to multiply u
vresult vector

References matvec().

Referenced by itersolve::itersolver(), and vankaprec().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_mergematrix ( type(matrix_t), pointer  A,
type(matrix_t), pointer  B,
integer, dimension(:), optional, pointer  PermA,
integer, dimension(:), optional, pointer  PermB 
)

Add another matrix B to matrix A, and eliminate B.

Parameters
aStructure holding the master matrix
bStructure holding the slave matrix
permaPermutation of the master dofs
permbPermutation of the slave dofs

References messages::fatal(), and messages::info().

Referenced by solverutils::mortarrobinsolver().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_moverow ( type(matrix_t)  A,
integer, intent(in)  n1,
integer, intent(in)  n2,
real(kind=dp), optional  coeff 
)

Add a row together with another row of a CRS matrix, and thereafter zero it.

Parameters
aStructure holding the matrix
[in]n1Row number to be copied and zerod
[in]n2Row number to be added
coeffOptional coefficient to multiply the row to be copied with

References crs_addtomatrixelement().

Referenced by solverutils::moverow().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_printmatrix ( type(matrix_t)  A)

Prints the values of the CRS matrix to standard output.

Parameters
aStructure holding matrix

Referenced by defutils::defaultsolve().

Here is the caller graph for this function:

subroutine crsmatrix::crs_printrhs ( type(matrix_t)  A)

Prints the values of the right-hand-side vector to standard output.

Parameters
aStructure holding matrix

Referenced by defutils::defaultsolve().

Here is the caller graph for this function:

subroutine crsmatrix::crs_removezeros ( type(matrix_t)  A)

Removes zeros from the matrix structure. This might be done in order to save memory, or to speed up the matrix operations. One must be carefull since the fact the an entry is zero does not always imply that it would be zero throughout the simulation.

Parameters
aThe matrix which will be returned with the non-zeros removed

References messages::info().

Here is the call graph for this function:

real(kind=dp) function crsmatrix::crs_rowsum ( type(matrix_t), intent(in)  A,
integer, intent(in)  k 
)

Computes the rowsoum of a given row in a CRS matrix.

Parameters
[in]aStructure holding matrix
[in]kRow index
subroutine crsmatrix::crs_rowsuminfo ( type(matrix_t), intent(in)  A,
real(kind=dp), dimension(:), optional, pointer  Values 
)

Computes information on the matrix rowsums.

References messages::info().

Referenced by crs_fctloworder().

Here is the call graph for this function:

Here is the caller graph for this function:

integer function crsmatrix::crs_search ( integer  N,
integer, dimension(:)  Array,
integer  VALUE 
)

Search CRS matrix for what?

Referenced by crs_addtomatrixelement(), crs_changematrixelement(), crs_getmatrixelement(), crs_setmatrixelement(), and crs_setsymmdirichlet().

Here is the caller graph for this function:

subroutine crsmatrix::crs_setmatrixelement ( type(matrix_t)  A,
integer, intent(in)  i,
integer, intent(in)  j,
real(kind=dp), intent(in)  VALUE 
)

Set a given value to an element of a CRS format matrix.

Parameters
aStructure holding the matrix
[in]irow number of the matrix element
[in]jcolumn number of the matrix element
[in]valuenew value of the matrix element

References crs_search().

Referenced by freesurface::poissonsolve(), and solverutils::setmatrixelement().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_setsymmdirichlet ( type(matrix_t)  A,
real(kind=dp), dimension(:)  b,
integer, intent(in)  n,
real(kind=dp), intent(in)  val 
)

When Dirichlet conditions are set by zeroing the row except for setting the diagonal entry to one, the matrix symmetry is broken. This routine maintains the symmetric structure of the matrix equation.

Parameters
aStructure holding matrix
[in]nIndex of the dofs to be fixed
bright-hand-side of the matrix equation
[in]valDirichlet value to be set

References crs_search(), and crs_zerorow().

Referenced by bulkassembly(), defutils::defaultdirichletbcs(), fetisolve::feti(), outletcompute(), solverutils::setdirichletpoint(), setdirichletpoints(), magnetodynamicsutils::setdoftovaluec(), magnetodynamicsutils::setdoftovaluer(), setelementvalues(), setlimitervalues(), and setpointvalues().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_sortbasicmatrix ( type(basicmatrix_t), target  A,
logical, intent(in), optional  ValuesToo 
)

Sort columns to ascending order for rows of a CRS format matrix- Optionally also sort the corresponging values of the matrix. This operates just on a basic matrix type.

Parameters
aStructure holding the matrix
[in]valuestooShould the values be sorted as well?

References sort(), and generalutils::sortf().

Referenced by sparitersolve::combinecrsmatindices().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_sortmatrix ( type(matrix_t)  A,
logical, intent(in), optional  ValuesToo 
)

Sort columns to ascending order for rows of a CRS format matrix- Optionally also sort the corresponging values of the matrix.

Parameters
aStructure holding the matrix
[in]valuestooShould the values be sorted as well?

References sort(), and generalutils::sortf().

Referenced by boundaryassembly(), bsolver(), elementutils::creatematrix(), crs_matrixmatrixmultiply(), eliminatedirichlet(), gebhardtfactors(), elementutils::initializematrix(), listmatrix::list_tocrs(), listmatrix::list_tocrsmatrix(), magnetodynamics2d(), radiationfactors(), rigidbody(), setperiodicboundariespass1(), and sparitersolve::sparinitsolve().

Here is the call graph for this function:

Here is the caller graph for this function:

type(matrix_t) function, pointer crsmatrix::crs_transpose ( type(matrix_t), pointer  A)

Calculate transpose of A in CRS format: B = A^T.

References generalutils::allocatematrix(), messages::error(), and messages::fatal().

Referenced by multigrid::amgsolve().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine crsmatrix::crs_transposematrixvectormultiply ( type(matrix_t), intent(in)  A,
real(kind=dp), dimension(*), intent(in)  u,
real(kind=dp), dimension(*), intent(out)  v 
)

Matrix vector product (v = A^T u) for a matrix given in CRS format.

Parameters
[in]uVector to be multiplied
[out]vResult vector
[in]aStructure holding matrix

Referenced by solverutils::transposematrixvectormultiply().

Here is the caller graph for this function:

subroutine crsmatrix::crs_zeromatrix ( type(matrix_t)  A)

Zero a CRS format matrix.

Parameters
aStructure holding the matrix

Referenced by acousticssolver(), gebhardtfactors(), solverutils::initializetozero(), freesurface::poissonsolve(), radiationfactors(), and stokessolver().

Here is the caller graph for this function:

subroutine crsmatrix::crs_zerorow ( type(matrix_t)  A,
integer  n 
)

Zero given row from a CRS format matrix.

Parameters
aStructure holding the matrix
nRow number to be zeroed

References lists::listgetlogical().

Referenced by crs_setsymmdirichlet(), freesurface::poissonsolve(), and solverutils::zerorow().

Here is the call graph for this function:

Here is the caller graph for this function:


The documentation for this module was generated from the following file: