ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_matrixMulTri Module Reference

This module contains the procedures for multiplication of a square triangular matrix in various transpositions with a general matrix.
More...

Data Types

interface  setMatMulTri
 Return the matrix solution to the system of linear equations with a triangular coefficient matrix \(\ms{C}\). More...
 

Variables

character(*, SK), parameter MODULE_NAME = "@pm_matrixMulTri"
 

Detailed Description

This module contains the procedures for multiplication of a square triangular matrix in various transpositions with a general matrix.

Such matrix multiplications frequently occur in solving general or triangular systems of equations with single or multiple right-hand sides (as a matrix).
Each system of linear equation is comprised of

  1. a coefficient matrix \(A\), which is assumed to be triangular in this module.
  2. the matrix of unknowns \(X\) for which the system must be solved,
  3. the intercept matrix \(B\) on the right-hand-side of the system of equations, which will be overwritten by the solution within this module.
Note
For general matrix-matrix or matrix-vector multiplications, see pm_matrixMulAdd.

The procedures under the generic interface setMatMulTri return the solution of the following corresponding linear triangular system of equations,

System of Equations Solution (corresponding matrix multiplication) BLAS routines
\(\ms{C}^{-1} \ms{X}^{\hphantom{-T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{C}^\hphantom{T} \ms{I} \) \(\ms{?trmv}~/~\ms{?trmm}\)
\(\ms{C}^{-T} \ms{X}^{\hphantom{-T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{C}^{T} \ms{I} \) \(\ms{?trmv}~/~\ms{?trmm}\)
\(\ms{C}^{-H} \ms{X}^{\hphantom{-T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{C}^{H} \ms{I} \) \(\ms{?trmv}~/~\ms{?trmm}\)
\(\ms{X}^{\hphantom{-T}} \ms{C}^{-1} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{I}^{\hphantom{T}} \ms{C} \) \(\ms{?trmv}~/~\ms{?trmm}\)
\(\ms{X}^{\hphantom{-T}} \ms{C}^{-T} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{I}^{\hphantom{T}} \ms{C}^{T} \) \(\ms{?trmv}~/~\ms{?trmm}\)
\(\ms{X}^{\hphantom{-T}} \ms{C}^{-H} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{I}^{\hphantom{T}} \ms{C}^{H} \) \(\ms{?trmv}~/~\ms{?trmm}\)
System of Equations Solution (corresponding matrix multiplication) BLAS routines
\(\ms{C}^{\hphantom{T}} \ms{X}^{\hphantom{T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{C}^{-1} \ms{I} \) \(\ms{?trsv}~/~\ms{?trsm}\)
\(\ms{C}^T \ms{X}^{\hphantom{T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{C}^{-T} \ms{I} \) \(\ms{?trsv}~/~\ms{?trsm}\)
\(\ms{C}^H \ms{X}^{\hphantom{T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{C}^{-H} \ms{I} \) \(\ms{?trsv}~/~\ms{?trsm}\)
\(\ms{X}^{\hphantom{T}} \ms{C}^{\hphantom{T}} = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{I}^{\hphantom{-T}} \ms{C}^{-1} \) \(\ms{?trsv}~/~\ms{?trsm}\)
\(\ms{X}^{\hphantom{T}} \ms{C}^T = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{I}^{\hphantom{-T}} \ms{C}^{-T} \) \(\ms{?trsv}~/~\ms{?trsm}\)
\(\ms{X}^{\hphantom{T}} \ms{C}^H = \alpha\ms{I}\) \( \ms{I} \leftarrow \alpha \ms{I}^{\hphantom{-T}} \ms{C}^{-H} \) \(\ms{?trsv}~/~\ms{?trsm}\)

where

  1. \(\cdot^-\) represents a matrix inversion.
  2. \(\cdot^T\) represents a Symmetric transpose.
  3. \(\cdot^H\) represents a Hermitian transpose.
  4. \(\ms{I}\) represents the Right-Hand-Side (RHS) of the system of linear equations.
  5. \(\ms{C}\) represents the triangular matrix of coefficients of the system of linear equations.
  6. \(\ms{X}\) represents the unknown matrix for which the system of linear equations must be solved (or equivalently, the result of the corresponding matrix multiplication).

Computational efficiency and practical considerations require the (matrix multiplication) solution to be written to \(\ms{I}\) instead of being returned explicitly to a separate argument as \(\ms{X}\).

BLAS/LAPACK equivalent:
The procedures under discussion combine, modernize, and extend the interface and functionalities of Version 3.11 of BLAS/LAPACK routine(s): STRMV, DTRMV, CTRMV, and ZTRMV, STRSV, DTRSV, CTRSV, and ZTRSV, STRMM, DTRMM, CTRMM, and ZTRMM, STRSM, DTRSM, CTRSM, and ZTRSM.
See also
nothing
inversion
lowerDiag
upperDiag
lowerUnit
upperUnit
transSymm
transHerm
transOrth
transUnit
setMatMulAdd
Todo:
High Priority: The BLAS routines ?TPMV, ?TPSV, ?TBMV and ?TBSV for Linear Full and Band packing format of the triangular input matrix must be implemented.
Todo:
Normal Priority: A benchmark comparison of the procedures of setMatMulTri with the optimized BLAS routine calls would be informative.
Test:
test_pm_matrixMul


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Friday 1:54 AM, April 21, 2017, Institute for Computational Engineering and Sciences (ICES), The University of Texas, Austin, TX

Variable Documentation

◆ MODULE_NAME

character(*, SK), parameter pm_matrixMulTri::MODULE_NAME = "@pm_matrixMulTri"

Definition at line 119 of file pm_matrixMulTri.F90.