SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
|
#include "superlu_sdefs.h"
Functions | |
void | psGetDiagU (int_t n, sLUstruct_t *LUstruct, gridinfo_t *grid, float *diagU) |
Copyright (c) 2003, The Regents of the University of California, through Lawrence Berkeley National Laboratory (subject to receipt of any required approvals from U.S. Dept. of Energy)
All rights reserved.
The source code is distributed under BSD license, see the file License.txt at the top-level directory.
void psGetDiagU | ( | int_t | n, |
sLUstruct_t * | LUstruct, | ||
gridinfo_t * | grid, | ||
float * | diagU | ||
) |
Purpose ======= GetDiagU extracts the main diagonal of matrix U of the LU factorization. Arguments ========= n (input) int Dimension of the matrix. LUstruct (input) sLUstruct_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 psgstrs_Bglobal.c.