SuperLU_DIST  4.0
superlu_dist on CPU and GPU clusters
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Macros | Functions
pdgstrs_Bglobal.c File Reference

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_ddefs.h"

Macros

#define ISEND_IRECV
 

Functions

void pdgstrs_Bglobal (int_t n, LUstruct_t *LUstruct, gridinfo_t *grid, double *B, int_t ldb, int nrhs, SuperLUStat_t *stat, int *info)
 

Detailed Description

Solves a system of distributed linear equations A*X = B with a general N-by-N matrix A using the LU factorization.

– 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

Macro Definition Documentation

#define ISEND_IRECV

Function Documentation

void pdgstrs_Bglobal ( int_t  n,
LUstruct_t LUstruct,
gridinfo_t grid,
double *  B,
int_t  ldb,
int  nrhs,
SuperLUStat_t stat,
int *  info 
)
Purpose
=======
pdgstrs_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 pdgstrf.
Arguments
=========
n      (input) int (global)
       The order of the system of linear equations.
LUstruct (input) LUstruct_t*
       The distributed data structures storing L and U factors.
       The L and U factors are obtained from pdgstrf for
       the possibly scaled and permuted matrix A.
       See superlu_ddefs.h for the definition of 'LUstruct_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) double*
       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