SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
psgstrs1.c File Reference

Solves a system of distributed linear equations. More...

#include "superlu_sdefs.h"
Include dependency graph for psgstrs1.c:

Macros

#define ISEND_IRECV
 

Functions

void psgstrs1 (superlu_dist_options_t *options, int_t n, sLUstruct_t *LUstruct, gridinfo_t *grid, float *x, int nrhs, SuperLUStat_t *stat, int *info)
 

Detailed Description

Solves a system of distributed linear equations.

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

Macro Definition Documentation

◆ ISEND_IRECV

#define ISEND_IRECV

Function Documentation

◆ psgstrs1()

void psgstrs1 ( superlu_dist_options_t options,
int_t  n,
sLUstruct_t LUstruct,
gridinfo_t grid,
float *  x,
int  nrhs,
SuperLUStat_t stat,
int *  info 
)
Purpose
=======

PSGSTRS1 solves a system of distributed linear equations

                  op( sub(A) ) * X = sub( B )

with a general N-by-N distributed matrix sub( A ) using the LU
factorization computed by PSGSTRF.

This routine is used only in the iterative refinement routine
psgsrfs_ABXglobal, assuming that the right-hand side is already
distributed in the diagonal processes.

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 to store L and U factors,
       and the permutation vectors.
       See superlu_sdefs.h for the definition of 'sLUstruct_t' structure.

grid   (input) gridinfo_t*
       The 2D process mesh.

x      (input/output) float*
       On entry, the right hand side matrix.
       On exit, the solution matrix if info = 0;

       NOTE: the right-hand side matrix is already distributed on
             the diagonal processes.

nrhs   (input) int (global)
       Number of right-hand sides.

stat   (output) SuperLUStat_t*
       Record the statistics about the triangular solves;
       See SuperLUStat_t structure defined in util.h.

info   (output) int*
       = 0: successful exit
    < 0: if info = -i, the i-th argument had an illegal value
Here is the call graph for this function: