SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
|
Solves a system of distributed linear equations A*X = B with a general N-by-N matrix A using the LU factorization. More...
#include "superlu_sdefs.h"
Macros | |
#define | ISEND_IRECV |
Functions | |
static void | gather_diag_to_all (int_t, int_t, float[], Glu_persist_t *, sLocalLU_t *, gridinfo_t *, int_t, int_t[], int_t[], float[], int_t, float[]) |
void | psgstrs_Bglobal (superlu_dist_options_t *options, int_t n, sLUstruct_t *LUstruct, gridinfo_t *grid, float *B, int_t ldb, int nrhs, SuperLUStat_t *stat, int *info) |
Solves a system of distributed linear equations A*X = B with a general N-by-N matrix A using the LU factorization.
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.
-- Distributed SuperLU routine (version 2.3) -- Lawrence Berkeley National Lab, Univ. of California Berkeley. October 15, 2008 Modified: Feburary 7, 2001 use MPI_Isend/MPI_Irecv October 2, 2001 use MPI_Isend/MPI_Irecv with MPI_Test October 15, 2008 use fewer MPI_Reduce
#define ISEND_IRECV |
|
static |
void psgstrs_Bglobal | ( | superlu_dist_options_t * | options, |
int_t | n, | ||
sLUstruct_t * | LUstruct, | ||
gridinfo_t * | grid, | ||
float * | B, | ||
int_t | ldb, | ||
int | nrhs, | ||
SuperLUStat_t * | stat, | ||
int * | info | ||
) |
Purpose ======= psgstrs_Bglobal solves a system of distributed linear equations A*X = B with a general N-by-N matrix A using the LU factorization computed by psgstrf. Arguments ========= options (input) superlu_dist_options_t* The structure defines the input parameters to control how the LU decomposition and triangular solve are performed. n (input) int (global) The order of the system of linear equations. LUstruct (input) sLUstruct_t* The distributed data structures storing L and U factors. The L and U factors are obtained from psgstrf for the possibly scaled and permuted matrix A. See superlu_ddefs.h for the definition of 'sLUstruct_t'. 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. Grid can be initialized by subroutine SUPERLU_GRIDINIT. See superlu_ddefs.h for the definition of 'gridinfo_t'. B (input/output) float* On entry, the right-hand side matrix of the possibly equilibrated and row permuted system. On exit, the solution matrix of the possibly equilibrated and row permuted system if info = 0; NOTE: Currently, the N-by-NRHS matrix B must reside on all processes when calling this routine. ldb (input) int (global) Leading dimension of matrix B. nrhs (input) int (global) Number of right-hand sides. stat (output) SuperLUStat_t* Record the statistics about the triangular solves. See util.h for the definition of 'SuperLUStat_t'. info (output) int* = 0: successful exit < 0: if info = -i, the i-th argument had an illegal value