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

Several matrix utilities. More...

#include <math.h>
#include "superlu_sdefs.h"
Include dependency graph for psutil.c:

Functions

int psCompRow_loc_to_CompCol_global (int_t need_value, SuperMatrix *A, gridinfo_t *grid, SuperMatrix *GA)
 Gather A from the distributed compressed row format to global A in compressed column format. More...
 
int psPermute_Dense_Matrix (int_t fst_row, int_t m_loc, int_t row_to_proc[], int_t perm[], float X[], int ldx, float B[], int ldb, int nrhs, gridinfo_t *grid)
 Permute the distributed dense matrix: B <= perm(X). perm[i] = j means the i-th row of X is in the j-th row of B. More...
 
void sLUstructInit (const int_t n, sLUstruct_t *LUstruct)
 Allocate storage in LUstruct. More...
 
void sLUstructFree (sLUstruct_t *LUstruct)
 Deallocate LUstruct. More...
 
void sDestroy_LU (int_t n, gridinfo_t *grid, sLUstruct_t *LUstruct)
 Destroy distributed L & U matrices. More...
 
int_t psgstrs_init (int_t n, int_t m_loc, int_t nrhs, int_t fst_row, int_t perm_r[], int_t perm_c[], gridinfo_t *grid, Glu_persist_t *Glu_persist, sSOLVEstruct_t *SOLVEstruct)
 
int sSolveInit (superlu_dist_options_t *options, SuperMatrix *A, int_t perm_r[], int_t perm_c[], int_t nrhs, sLUstruct_t *LUstruct, gridinfo_t *grid, sSOLVEstruct_t *SOLVEstruct)
 Initialize the data structure for the solution phase. More...
 
void sSolveFinalize (superlu_dist_options_t *options, sSOLVEstruct_t *SOLVEstruct)
 Release the resources used for the solution phase. More...
 
void sDestroy_A3d_gathered_on_2d (sSOLVEstruct_t *SOLVEstruct, gridinfo3d_t *grid3d)
 
void psinf_norm_error (int iam, int_t n, int_t nrhs, float x[], int_t ldx, float xtrue[], int_t ldxtrue, MPI_Comm slucomm)
 Check the inf-norm of the error vector. More...
 
void sDestroy_Tree (int_t n, gridinfo_t *grid, sLUstruct_t *LUstruct)
 Destroy broadcast and reduction trees used in triangular solve. More...
 

Detailed Description

Several matrix utilities.

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.0) --
Lawrence Berkeley National Lab, Univ. of California Berkeley.
March 15, 2003

Last modified:
    December 28, 2022

Function Documentation

◆ psCompRow_loc_to_CompCol_global()

int psCompRow_loc_to_CompCol_global ( int_t  need_value,
SuperMatrix A,
gridinfo_t grid,
SuperMatrix GA 
)

Gather A from the distributed compressed row format to global A in compressed column format.

Here is the call graph for this function:

◆ psgstrs_init()

int_t psgstrs_init ( int_t  n,
int_t  m_loc,
int_t  nrhs,
int_t  fst_row,
int_t  perm_r[],
int_t  perm_c[],
gridinfo_t grid,
Glu_persist_t Glu_persist,
sSOLVEstruct_t SOLVEstruct 
)
Purpose
=======
  Set up the communication pattern for redistribution between B and X
  in the triangular solution.

Arguments
=========

n      (input) int (global)
       The dimension of the linear system.

m_loc  (input) int (local)
       The local row dimension of the distributed input matrix.

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

fst_row (input) int (global)
       The row number of matrix B's first row in the global matrix.

perm_r (input) int* (global)
       The row permutation vector.

perm_c (input) int* (global)
       The column permutation vector.

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

◆ psinf_norm_error()

void psinf_norm_error ( int  iam,
int_t  n,
int_t  nrhs,
float  x[],
int_t  ldx,
float  xtrue[],
int_t  ldxtrue,
MPI_Comm  slucomm 
)

Check the inf-norm of the error vector.

◆ psPermute_Dense_Matrix()

int psPermute_Dense_Matrix ( int_t  fst_row,
int_t  m_loc,
int_t  row_to_proc[],
int_t  perm[],
float  X[],
int  ldx,
float  B[],
int  ldb,
int  nrhs,
gridinfo_t grid 
)

Permute the distributed dense matrix: B <= perm(X). perm[i] = j means the i-th row of X is in the j-th row of B.

Here is the call graph for this function:

◆ sDestroy_A3d_gathered_on_2d()

void sDestroy_A3d_gathered_on_2d ( sSOLVEstruct_t SOLVEstruct,
gridinfo3d_t grid3d 
)

◆ sDestroy_LU()

void sDestroy_LU ( int_t  n,
gridinfo_t grid,
sLUstruct_t LUstruct 
)

Destroy distributed L & U matrices.

Here is the call graph for this function:

◆ sDestroy_Tree()

void sDestroy_Tree ( int_t  n,
gridinfo_t grid,
sLUstruct_t LUstruct 
)

Destroy broadcast and reduction trees used in triangular solve.

Here is the call graph for this function:

◆ sLUstructFree()

void sLUstructFree ( sLUstruct_t LUstruct)

Deallocate LUstruct.

◆ sLUstructInit()

void sLUstructInit ( const int_t  n,
sLUstruct_t LUstruct 
)

Allocate storage in LUstruct.

Here is the call graph for this function:

◆ sSolveFinalize()

void sSolveFinalize ( superlu_dist_options_t options,
sSOLVEstruct_t SOLVEstruct 
)

Release the resources used for the solution phase.

Here is the call graph for this function:

◆ sSolveInit()

int sSolveInit ( superlu_dist_options_t options,
SuperMatrix A,
int_t  perm_r[],
int_t  perm_c[],
int_t  nrhs,
sLUstruct_t LUstruct,
gridinfo_t grid,
sSOLVEstruct_t SOLVEstruct 
)

Initialize the data structure for the solution phase.

Here is the call graph for this function: