You cannot select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
psblas3/docs/pdf/psbrout.tex

870 lines
28 KiB
TeX

19 years ago
\section{Algebraic routines}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% DENSE MATRIX SUM
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_geaxpby}{General Dense Matrix Sum}
19 years ago
This subroutine is an interface to the computational kernel for
dense matrix sum:
\[ y \leftarrow \alpha\> x+ \beta y \]
where:
\begin{description}
\item[$x$] represents the global dense submatrix $x_{:, jx:jx+n-1}$
\item[$y$] represents the global dense submatrix $y_{:, jy:jy+n-1}$
\end{description}
\syntax{call psb\_geaxpby}{alpha, x, beta, y, desc\_a, info}
\syntax*{call psb\_geaxpby}{alpha, x, beta, y, desc\_a, info, n, jx, jy}
19 years ago
%( calculating y <- alpha*x+beta*y )
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$x$, $y$, $\alpha$, $\beta$ & {\bf Subroutine}\\
\hline
Single Precision Real & psb\_geaxpby \\
Long Precision Real & psb\_geaxpby \\
Long Precision Complex & psb\_geaxpby \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90axpby}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[alpha] the scalar $\alpha$.\\
Scope: {\bf global} \\
Type: {\bf required} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90axpby}.
\item[x] the local portion of global dense matrix
$x$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type
specified in Table~\ref{tab:f90axpby}. The rank of $x$ must be the same of $y$.
\item[beta] the scalar $\beta$.\\
Scope: {\bf global} \\
Type: {\bf required} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90axpby}.
\item[y] the local portion of the global dense matrix
$y$. \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER
attributecontaining numbers of the type
indicated in Table~\ref{tab:f90axpby}. The rank of $y$ must be the same of $x$.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[n] number of columns in dense submatrices $x$ and $y$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ and $y$ are of rank 2.\\
Default: \verb|min(size(x,2),size(y,2))|.\\
Specified as: an integer variable $n\ge 0$.
\item[jx] the column index of the global dense matrix $x$,
identifying the first column of the submatrix $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ and $y$ are of rank 2.\\
Default: $jx = 1$.\\
Specified as: an integer variable $jx\ge 1$.
\item[jy] the column index of the global dense matrix $y$,
identifying the first column of the submatrix $y$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ and $y$ are of rank 2.\\
Default: $jy = 1$.\\
Specified as: an integer variable $jy\ge 1$.
\end{description}
\begin{description}
\item[\bf On Return]
\item[y] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array containing numbers of the type
indicated in Table~\ref{tab:f90axpby}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% F90DOT PRODUCT
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_gedot}{Dot Product}
19 years ago
This function computes dot product between two vectors $x$ and
$y$.\\
If $x$ and $y$ are double precision real or single precision real vectors
computes dot-product as:
\[dot \leftarrow x^T y\]
Else if $x$ and $y$ are double precision complex vectors then computes dot-product as:
\[dot \leftarrow x^H y\]
where:
\begin{description}
\item[$x$] represents the global subvector $x_{:,jx}$
\item[$y$] represents the global subvector $y_{:,jy}$
\end{description}
\syntax{psb\_gedot}{x, y, desc\_a, info}
\syntax*{psb\_gedot}{x, y, desc\_a, info, jx, jy}
19 years ago
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$dot$, $x$, $y$ & {\bf Function}\\
\hline
Single Precision Real & psb\_gedot\\
Long Precision Real & psb\_gedot \\
Long Precision Complex & psb\_gedot \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90dot}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[x] the local portion of global dense matrix
$x$. This function computes the location of the first element of
local subarray used, based on $jx$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a pointer to array of rank one or two
containing numbers of type specified in
Table~\ref{tab:f90dot}. The rank of $x$ must be the same of $y$.
\item[y] the local portion of global dense matrix
$y$. This function computes the location of the first element of
local subarray used, based on $iy, jy$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a pointer to array of rank one or two
containing numbers of type specified in
Table~\ref{tab:f90dot}. The rank of $y$ must be the same of $x$.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[jx] the column index of global dense matrix $x$,
identifying the column of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ and $y$ are of rank 2.\\
Default: $jx = 1$.\\
\item[jy] the column index of global dense matrix $y$,
identifying the column of subvector $y$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ and $y$ are of rank 2.\\
Default: $jy = 1$.\\
Specified as: an integer variable $jy\ge 1$.
\item[\bf On Return]
\item[Function value] is the dot product of subvectors $x$ and $y$.\\
Scope: {\bf global} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90dot}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% F90DOT PRODUCT
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_gedot}{Generalized Dot Product}
19 years ago
This subroutine computes a series of dot products among the columns of
two dense matrices $x$ and $y$:
\[ res(i) \leftarrow x(:,i)^T y(:,i)\]
If the matrices are complex, then the
usual convention applies, i.e. the conjugate transpose of $x$ is
used. If $x$ and $y$ are of rank one, then $res$ is a scalar, else it
is a rank one array.
\syntax{psb\_gedot}{res, x, y, desc\_a, info}
19 years ago
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$res$, $x$, $y$ & {\bf Subroutine}\\
\hline
Single Precision Real & psb\_gedot\\
Long Precision Real & psb\_gedot \\
Long Precision Complex & psb\_gedot \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90mdot}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[x] the local portion of global dense matrix
$x$. \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a pointer to array of rank one or two
containing numbers of type specified in
Table~\ref{tab:f90mdot}. The rank of $x$ must be the same of $y$.
\item[y] the local portion of global dense matrix
$y$. \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a pointer to array of rank one or two
containing numbers of type specified in
Table~\ref{tab:f90mdot}. The rank of $y$ must be the same of $x$.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[\bf On Return]
\item[res] is the dot product of subvectors $x$ and $y$.\\
Scope: {\bf global} \\
Specified as: a number or a rank-one array of the data type indicated
in Table~\ref{tab:f90dot}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% VECTOR INFINITY-NORM
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_geamax}{Infinity-Norm of Vector}
19 years ago
This function computes
the infinity-norm of a vector $x$.\\
If $x$ is double precision real or single precision real vector
computes infinity norm as:
\[ amax \leftarrow \max_i |x_i|\]
else if $x$ is double precision complex vector then computes infinity-norm as:
\[ amax \leftarrow \max_i {(|re(x_i)| + |im(x_i)|)}\]
where:
\begin{description}
\item[$x$] represents the global subvector $x_{:,jx}$
\end{description}
\syntax{psb\_geamax}{x, desc\_a, info}
\syntax*{psb\_geamax}{x, desc\_a, info, jx}
19 years ago
\begin{table}[h]
\begin{center}
\begin{tabular}{lll}
\hline
$amax$ & $x$ & {\bf Function}\\
\hline
Single Precision Real&Single Precision Real & psb\_geamax\\
Long Precision Real&Long Precision Real & psb\_geamax \\
Long Precision Real&Long Precision Complex & psb\_geamax \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90amax}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[x] the local portion of global dense matrix
$x$. This function computes the location of the first element of
local subarray used, based on $jx$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90amax}.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[jx] the column index of global dense matrix $x$,
identifying the column of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ is of rank 2.\\
Default: $jx = 1$\\
Specified as: an integer variable $jx\ge 1$.
\item[\bf On Return]
\item[Function value] is the infinity norm of subvector $x$.\\
Scope: {\bf global} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90amax}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% Infinity norm
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_geamax}{Generalized Infinity Norm}
19 years ago
This subroutine computes a series of infinity norms on the columns of
a dense matrix $x$:
\[ res(i) \leftarrow \max_k |x(k,i)| \]
\syntax{psb\_geamax}{res, x, desc\_a, info}
19 years ago
\begin{table}[h]
\begin{center}
\begin{tabular}{lll}
\hline
$res$& $x$& {\bf Subroutine}\\
\hline
Single Precision Real &Single Precision Real & psb\_geamax\\
Long Precision Real &Long Precision Real & psb\_geamax\\
Long Precision Real &Long Precision Complex & psb\_geamax\\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90mamax}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[x] the local portion of global dense matrix
$x$. \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90mamax}.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[\bf On Return]
\item[res] is the infinity norm of the columns of $x$.\\
Scope: {\bf global} \\
Specified as: a number or a rank-one array of the data type indicated
in Table~\ref{tab:f90amax}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% 1-NORM OF A VECTOR
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_geasum}{1-Norm of Vector}
19 years ago
This function computes the 1-norm of a vector $x$.\\
If $x$ is double precision real or single precision real vector
computes 1-norm as:
\[ asum \leftarrow \|x_i\|\]
else if $x$ ic double precision complex vector then computes 1-norm as:
\[ asum \leftarrow \|re(x)\|_1 + \|im(x)\|_1\]
where:
\begin{description}
\item[$x$] represents the global subvector $x_{:,jx}$
\end{description}
\syntax{psb\_geasum}{x, desc\_a, info}
\syntax*{psb\_geasum}{x, desc\_a, info, jx}
19 years ago
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$dot$, $x$, $y$ & {\bf Function}\\
\hline
Single Precision Real & psb\_geasum\\
Long Precision Real & psb\_geasum \\
Long Precision Complex & psb\_geasum \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90asum}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[x] the local portion of global dense matrix
$x$. This function computes the location of the first element of
local subarray used, based on $jx$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90asum}.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[jx] the column index of global dense matrix $x$,
identifying the column of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ is of rank 2.\\
Default: $jx = 1$\\
Specified as: an integer variable $jx\ge 1$.
\item[\bf On Return]
\item[Function value] is the 1-norm of subvector $x$.\\
Scope: {\bf global} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90asum}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% 2-NORM OF A VECTOR
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine {psb\_genrm2}{2-Norm of Vector}
19 years ago
This function computes the 2-norm of a vector $x$.\\
If $x$ is double precision real or single precision real vector
computes 2-norm as:
\[ nrm2 \leftarrow \sqrt{x^T x}\]
else if $x$ is double precision complex vector then computes 2-norm as:
\[ nrm2 \leftarrow \sqrt{x^H x}\]
where:
\begin{description}
\item[$x$] represents the global subvector $x_{:,jx}$
\end{description}
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$nrm2$, $x$ & {\bf Function}\\
\hline
Single Precision Real & psb\_genrm2\\
Long Precision Real & psb\_genrm2 \\
Long Precision Complex & psb\_genrm2 \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90nrm2}}
\end{table}
\syntax{psb\_genrm2}{x, desc\_a, info}
\syntax*{psb\_genrm2}{x, desc\_a, info, jx}
19 years ago
\begin{description}
\item[\bf On Entry]
\item[x] the local portion of global dense matrix
$x$. This function computes the location of the first element of
local subarray used, based on $jx$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90nrm2}.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[jx] the column index of global dense matrix $x$,
identifying the column of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ is of rank 2.\\
Default: $jx = 1$\\
Specified as: an integer variable $jx\ge 1$.
\item[\bf On Return]
\item[Function Value] is the 2-norm of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf required} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90nrm2}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% INFINITY-NORM OF A MATRIX
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_spnrmi}{Infinity Norm of Sparse Matrix}
19 years ago
This function computes the infinity-norm of a matrix $A$:\\
\[ nrmi \leftarrow \|A\|_\infty \]
where:
\begin{description}
\item[$A$] represents the global matrix $A$
\end{description}
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$nrmi$, $A$ & {\bf Function}\\
\hline
Single Precision Real & psb\_spnrmi\\
Long Precision Real & psb\_spnrmi \\
Long Precision Complex & psb\_spnrmi \\
19 years ago
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90nrmi}}
\end{table}
\syntax{psb\_spnrmi}{A, desc\_a, info}
19 years ago
\begin{description}
\item[\bf On Entry]
\item[a] the local portion of the global sparse matrix
$A$. \\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \spdata.
19 years ago
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[\bf On Return]
\item[Function value] is the infinity-norm of sparse submatrix $A$.\\
Scope: {\bf global} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90nrmi}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% SPARSE MATRIX by DENSE MATRIX PRODUCT
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_spmm}{Sparse Matrix by Dense Matrix Product}
This subroutine computes the Sparse Matrix by Dense Matrix Product:
\begin{equation}
y \leftarrow \alpha P_r A P_c x + \beta y
\label{eq:f90spmm_no_tra}
\end{equation}
\begin{equation}
y \leftarrow \alpha P_r A^T P_c x + \beta y
\label{eq:f90spmm_tra}
\end{equation}
\begin{equation}
y \leftarrow \alpha P_r A^H P_c x + \beta y
\label{eq:f90spmm_con}
\end{equation}
where:
\begin{description}
\item[$x$] is the global dense submatrix $x_{:, jx:jx+k-1}$
\item[$y$] is the global dense submatrix $y_{:, jy:jy+k-1}$
\item[$A$] is the global sparse submatrix $A$
\item[$P_r, P_c$] are the permutation matrices.
\end{description}
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$A$, $x$, $y$, $\alpha$, $\beta$ & {\bf Subroutine}\\
\hline
Single Precision Real & psb\_spmm\\
Long Precision Real & psb\_spmm \\
Long Precision Complex & psb\_spmm \\
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90spmm}}
\end{table}
\syntax{CALL psb\_spmm}{alpha, a, x, beta, y, desc\_a, info}
\syntax*{CALL psb\_spmm}{alpha, a, x, beta, y,desc\_a, info,
trans, k, jx, jy, work}
\begin{description}
\item[\bf On Entry]
\item[alpha] the scalar $\alpha$.\\
Scope: {\bf global} \\
Type: {\bf required}\\
Specified as: a number of the data type indicated in
Table~\ref{tab:f90spmm}.
\item[a] the local portion of the sparse matrix
$A$. \\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \spdata.
19 years ago
\item[x] the local portion of global dense matrix
$x$. This subroutine computes the location of the first element of
local subarray used, based on $jx$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90spmm}. The rank of $x$ must be the same of $y$.
\item[beta] the scalar $\beta$.\\
Scope: {\bf global} \\
Type: {\bf required} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90spmm}.
\item[y] the local portion of global dense matrix
$y$. This subroutine computes the location of the first element of
local subarray used, based on $jy$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90spmm}. The rank of $y$ must be the same of $x$.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[trans] indicate what kind of operation to perform.
\begin{description}
\item[trans = N] the operation is specified by equation \ref{eq:f90spmm_no_tra}
\item[trans = T] the operation is specified by equation
\ref{eq:f90spmm_tra}
\item[trans = C] the operation is specified by equation
\ref{eq:f90spmm_con}
\end{description}
Scope: {\bf global} \\
Type: {\bf optional}\\
Default: $trans = N$\\
Specified as: a character variable.
\item[k] number of columns in dense submatrices $x$ and $y$. \\
Scope: {\bf global} \\
Type: {\bf optional}\\
Default: \verb|min(size(x,2)-jx+1,size(y,2)-jy+1)|\\
Specified as: an integer variable $ k \ge 1$.
\item[jx] the column index of global dense matrix $x$,
identifying the column of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ is of rank 2.\\
Default: $iy = 1$\\
Specified as: an integer variable $jx\ge 1$.
\item[jy] the column index of global dense matrix $y$,
identifying the column of subvector $y$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $y$ is of rank 2.\\
Default: $jy = 1$\\
Specified as: an integer variable $jy\ge 1$.
\item[work] the work array.\\
Scope: {\bf local} \\
Specified as: a rank one array of the same type of $x$ and $y$ with
the POINTER attribute.
\item[\bf On Return]
\item[y] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a pointer to array of rank one or two
containing numbers of type specified in
Table~\ref{tab:f90spmm}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%
% TRIANGULAR SYSTEM SOLVE
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subroutine{psb\_spsm}{Triangular System Solve}
This subroutine computes the Triangular System Solve:
\begin{eqnarray*}
y &\leftarrow& \alpha P_r T^{-1} P_c x + \beta y\\
y &\leftarrow& \alpha D P_r T^{-1} P_c x + \beta y\\
y &\leftarrow& \alpha P_r T^{-1} P_c D x + \beta y\\
y &\leftarrow& \alpha P_r T^{-T} P_c x + \beta y\\
y &\leftarrow& \alpha D P_r T^{-T} P_c x + \beta y\\
y &\leftarrow& \alpha P_r T^{-T} P_c D x + \beta y\\
y &\leftarrow& \alpha P_r T^{-H} P_c x + \beta y\\
y &\leftarrow& \alpha D P_r T^{-H} P_c x + \beta y\\
y &\leftarrow& \alpha P_r T^{-H} P_c D x + \beta y\\
\end{eqnarray*}
where:
\begin{description}
\item[$x$] is the global dense submatrix $x_{:, jx:jx+n-1}$
\item[$y$] is the global dense submatrix $y_{:, jy:jy+n-1}$
\item[$T$] is the global sparse block triangular submatrix $T$
\item[$D$] is the scaling diagonal matrix.
\item[$P_r, P_c$] are the permutation matrices.
\end{description}
\syntax{CALL psb\_spsm}{alpha, t, x, beta, y, desc\_a, info}
\syntax*{CALL psb\_spsm}{alpha, t, x, beta, y, desc\_a, info,
trans, unit, choice, diag, n, jx, jy, work}
\begin{table}[h]
\begin{center}
\begin{tabular}{ll}
\hline
$T$, $x$, $y$, $D$, $\alpha$, $\beta$ & {\bf Subroutine}\\
\hline
Single Precision Real & psb\_spsm\\
Long Precision Real & psb\_spsm \\
Long Precision Complex & psb\_spsm \\
\hline
\end{tabular}
\end{center}
\caption{Data types\label{tab:f90spsm}}
\end{table}
\begin{description}
\item[\bf On Entry]
\item[alpha] the scalar $\alpha$.\\
Scope: {\bf global} \\
Type: {\bf required}\\
Specified as: a number of the data type indicated in
Table~\ref{tab:f90spsm}.
\item[t] the global portion of the sparse matrix
$T$. \\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data type specified in
\S~\ref{sec:datastruct}.
\item[x] the local portion of global dense matrix
$x$. This subroutine computes the location of the first element of
local subarray used, based on $jx$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90spsm}. The rank of $x$ must be the same of $y$.
\item[beta] the scalar $\beta$.\\
Scope: {\bf global} \\
Type: {\bf required} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90spsm}.
\item[y] the local portion of global dense matrix
$y$. This subroutine computes the location of the first element of
local subarray used, based on $jy$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a rank one or two array with the POINTER attribute
containing numbers of type specified in
Table~\ref{tab:f90spsm}. The rank of $y$ must be the same of $x$.
\item[desc\_a] contains data structures for communications.\\
Scope: {\bf local} \\
Type: {\bf required}\\
Specified as: a structured data of type \descdata.
19 years ago
\item[trans] specify with {\em unitd} the operation to perform.
\begin{description}
\item[trans = 'N'] the operation is with no transposed matrix
\item[trans = 'T'] the operation is with transposed matrix.
\item[trans = 'C'] the operation is with conjugate transposed matrix.
\end{description}
Scope: {\bf global} \\
Type: {\bf optional}\\
Default: $trans = N$\\
Specified as: a character variable.
\item[unitd] specify with {\em trans} the operation to perform.
\begin{description}
\item[unitd = 'U'] the operation is with no scaling
\item[unitd = 'L'] the operation is with left scaling
\item[unitd = 'R'] the operation is with right scaling.
\end{description}
Scope: {\bf global} \\
Type: {\bf optional}\\
Default: $unitd = U$\\
Specified as: a character variable.
\item[choice] specify whether a cleanup of the overlapped elements is
required on exit.
\begin{description}
\item[choice = .false.] no cleanup on exit
\item[choice = .true.] cleanup on exit.\\
\end{description}
Scope: {\bf global} \\
Type: {\bf optional}\\
Default: $choice = .true.$\\
Specified as: a logical variable.
\item[diag] the diagonal scaling matrix.\\
Scope: {\bf local} \\
Type: {\bf optional}\\
Default: $diag(1) = 1 (no scaling)$\\
Specified as: a rank one array containing numbers of the type
indicated in Table~\ref{tab:f90spsm}.
\item[n] number of columns in dense submatrices $x$ and $y$. \\
Scope: {\bf global} \\
Type: {\bf optional}\\
Default: \verb|min(size(x,2)-jx+1,size(y,2)-jy+1)|\\
Specified as: an integer variable $ n \ge 0$.
\item[jx] the column index of global dense matrix $x$,
identifying the column of subvector $x$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $x$ is of rank 2.\\
Default: $jx = 1 $\\
Specified as: an integer variable $jx\ge 1$.
\item[jy] the column index of global dense matrix $y$,
identifying the column of subvector $y$.\\
Scope: {\bf global} \\
Type: {\bf optional}; can only be present if $y$ is of rank 2.\\
Default: $jy = 1 $\\
Specified as: an integer variable $jy\ge 1$. \\
Scope: {\bf global} \\
Specified as: a number of the data type indicated in Table~\ref{tab:f90spsm}.
\item[work] the work array. \\
Scope: {\bf local} \\
Type: {\bf optional}\\
Specified as: a rank one array of the same type of $x$ with the
POINTER attribute.
\item[\bf On Return]
\item[y] the local portion of global dense matrix
$y$. This subroutine computes the location of the first element of
local subarray used, based on $jy$ and the field $matrix\_data$ of $desc\_a$ . \\
Scope: {\bf local} \\
Type: {\bf required} \\
Specified as: a pointer to array of rank one or two
containing numbers of type specified in
Table~\ref{tab:f90spsm}.
\item[info] the local portion of result submatrix $y$.\\
Scope: {\bf local} \\
Type: {\bf required} \\
An integer value that contains an error code.
\end{description}
%%% Local Variables:
%%% mode: latex
%%% TeX-master: "userguide"
%%% End: