SuperLU Distributed 9.0.0
gpu3d
pddistribute3d.h
Go to the documentation of this file.
1#pragma once
2
3#include "superlu_ddefs.h"
4
5// Defining DIST3D temoprarily to track changes in this file while fully
6// integrating the 3D code
7// #define DIST3D
9 dScalePermstruct_t *ScalePermstruct,
10 Glu_freeable_t *Glu_freeable, dLUstruct_t *LUstruct,
11 gridinfo3d_t *grid3d);
12
13
15 Glu_freeable_t *Glu_freeable, int_t *xsup, int_t *supno,
16 gridinfo3d_t *grid3d, int_t *colptr[], int_t *rowind[],
17 double *a[]);
18
20 dLUstruct_t *LUstruct,
21 int_t *xa,
22 int_t *asub,
23 double *a,
25 gridinfo3d_t *grid3d,
26 int_t nsupers,
27 float *mem_use);
28int_t ComputeLDAspa_Ilsum( int_t nsupers, int_t* ilsum, gridinfo3d_t* grid3d) ;
29
30// void propagateDataThroughMatrixBlocks(int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t* grid3d, int_t *xusub, int_t *usub, int_t **ToSendR,
31// int_t *ToSendD, int_t *Urb_length, int_t *rb_marker, int_t *Urb_fstnz, int_t *Ucbs, int_t *ToRecv);
32
33void propagateDataThroughMatrixBlocks(int_t nsupers, Glu_freeable_t *Glu_freeable, dLUstruct_t *LUstruct, gridinfo3d_t* grid3d,
34int_t *Urb_length, int_t *rb_marker, int_t *Urb_fstnz, int_t *Ucbs,
35int **ToSendR, int *ToSendD, int *ToRecv);
37
39 dScalePermstruct_t *ScalePermstruct,
40 Glu_freeable_t *Glu_freeable,
41 dLUstruct_t *LUstruct, gridinfo3d_t *grid3d);
42
43void dnewTrfPartitionInit(int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d);
44
45
46int compareInt_t(void *a, void *b);
47int compareInt(void *a, void *b);
48int compareDouble(void *a, void *b);
49int dist_checkArrayEq(void *arr, int length, MPI_Datatype datatype, int src_rank, int dest_rank, MPI_Comm communicator, int (*compare)(void *, void *));
int_t checkDist3DLUStruct(dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
Definition: pddistribute-aux3d.c:351
int compareDouble(void *a, void *b)
Compares two doubles for equality.
Definition: distCheckArray.c:63
int compareInt(void *a, void *b)
Compares two integers for equality.
Definition: distCheckArray.c:42
int dist_checkArrayEq(void *arr, int length, MPI_Datatype datatype, int src_rank, int dest_rank, MPI_Comm communicator, int(*compare)(void *, void *))
Checks whether arrays at two MPI ranks are identical.
Definition: distCheckArray.c:93
float pddistribute3d(superlu_dist_options_t *options, int_t n, SuperMatrix *A, dScalePermstruct_t *ScalePermstruct, Glu_freeable_t *Glu_freeable, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t ComputeLDAspa_Ilsum(int_t nsupers, int_t *ilsum, gridinfo3d_t *grid3d)
int_t dReDistribute_A3d(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct, Glu_freeable_t *Glu_freeable, int_t *xsup, int_t *supno, gridinfo3d_t *grid3d, int_t *colptr[], int_t *rowind[], double *a[])
void propagateDataThroughMatrixBlocks(int_t nsupers, Glu_freeable_t *Glu_freeable, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, int_t *Urb_length, int_t *rb_marker, int_t *Urb_fstnz, int_t *Ucbs, int **ToSendR, int *ToSendD, int *ToRecv)
Definition: pddistribute-aux3d.c:281
int compareInt_t(void *a, void *b)
Compares two integers for equality.
Definition: distCheckArray.c:20
void dnewTrfPartitionInit(int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
Definition: d3DPartition.c:5
void dbcastPermutedSparseA(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct, Glu_freeable_t *Glu_freeable, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
Definition: d3DPartition.c:187
void dpropagate_A_to_LU3d(dLUstruct_t *LUstruct, int_t *xa, int_t *asub, double *a, superlu_dist_options_t *options, gridinfo3d_t *grid3d, int_t nsupers, float *mem_use)
Definition: pddistribute-aux3d.c:19
Definition: superlu_defs.h:506
Definition: supermatrix.h:54
Definition: superlu_ddefs.h:340
Definition: superlu_ddefs.h:76
Definition: superlu_defs.h:414
Definition: superlu_defs.h:728
Distributed SuperLU data types and function prototypes.
int64_t int_t
Definition: superlu_defs.h:119