#include "superlu_ddefs.h"
Functions/Subroutines | |
| void | GetDiagU (int_t n, LUstruct_t *LUstruct, gridinfo_t *grid, double *diagU) |
-- Auxiliary routine in distributed SuperLU (version 1.0) -- Lawrence Berkeley National Lab, Univ. of California Berkeley. Xiaoye S. Li April 16, 2002
| void GetDiagU | ( | int_t | n, | |
| LUstruct_t * | LUstruct, | |||
| gridinfo_t * | grid, | |||
| double * | diagU | |||
| ) |
Purpose =======
GetDiagU extracts the main diagonal of matrix U of the LU factorization.
Arguments =========
n (input) int
Dimension of the matrix.
LUstruct (input) LUstruct_t*
The data structures to store the distributed L and U factors.
see superlu_ddefs.h for its definition.
grid (input) gridinfo_t*
The 2D process mesh. It contains the MPI communicator, the number
of process rows (NPROW), the number of process columns (NPCOL),
and my process rank. It is an input argument to all the
parallel routines.
diagU (output) double*, dimension (n)
The main diagonal of matrix U.
On exit, it is available on all processes.
Note ====
The diagonal blocks of the L and U matrices are stored in the L data structures, and are on the diagonal processes of the 2D process grid.
This routine is modified from gather_diag_to_all() in pdgstrs_Bglobal.c.
1.5.5