HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
tridiagLU.h
Go to the documentation of this file.
1 
6 #ifndef _TRIDIAGLU_H_
7 #define _TRIDIAGLU_H_
8 
9 /*
10 
11  Parallel direct solver for tridiagonal systems
12 
13  tridiagLU (a,b,c,x,n,ns,r,m) - Parallel tridiagonal solver
14 
15  Solves the tridiagonal system in parallel by reordering the
16  points such that the first point of each subdomain is placed
17  at the end.
18 
19  The interior points are eliminated in parallel, resulting in
20  a reduced system consisting of the first point of each sub-
21  domain.
22 
23  This reduced system is solved either by the gather-and-
24  solve (tridiagLUGS) or the recursive-doubling (tridiagLURD)
25  algorithms.
26 
27  tridiagLUGS(a,b,c,x,n,ns,r,m) - Tridiagonal solver based on
28  "gather and solve"
29 
30  Each of the "ns" systems is gathered on one processor,
31  solved in serial, and the solution scattered back. The
32  parallelism lies in solving the "ns" different systems
33  by multiple processors (i.e., each processor solves
34  ~ns/nproc number of systems in serial).
35 
36  Arguments:-
37  a [0,ns-1]x[0,n-1] double* subdiagonal entries
38  b [0,ns-1]x[0,n-1] double* diagonal entries
39  c [0,ns-1]x[0,n-1] double* superdiagonal entries
40  x [0,ns-1]x[0,n-1] double* right-hand side (solution)
41  n int local size of the system
42  ns int number of systems to solve
43  r TridiagLU* structure containing paramters
44  for the tridiagonal solve and
45  the walltimes:
46  total_time
47  stage1_time
48  stage2_time
49  stage3_time
50  stage4_time
51  ** Note that these are process-specific. Calling
52  function needs to do something to add/average
53  them to get some global value.
54  ** Can be NULL if runtimes are not needed.
55  m MPI_Comm* MPI Communicator
56 
57  For a,b,c, and x, [0,ns-1] is the inner loop, i.e., the i-th row of the
58  d-th system is a[i*ns+d], b[i*ns+d], c[i*ns+d] and x[i*ns+d].
59 
60  Return value (int) -> 0 (successful solve), -1 (singular system)
61 
62  Note:-
63  x contains the final solution (right-hand side is replaced)
64  a,b,c are not preserved
65  On rank=0, a[0*ns+d] has to be zero for all d.
66  On rank=nproc-1, c[(n-1)*ns+d] has to be zero for all d.
67 
68  For a serial tridiagonal solver, compile with the flag "-Dserial"
69  or send NULL as the argument for the MPI communicator.
70 
71 */
72 
74 #define _TRIDIAG_JACOBI_ "jacobi"
75 
76 #define _TRIDIAG_GS_ "gather-and-solve"
77 
84 typedef struct _tridiagLU_ {
85 
86  /* Parameters for tridiagLU() */
87 
91  char reducedsolvetype[50];
92 
94  int maxiter;
95  double atol,
96  rtol;
97  int exititer;
98  double exitnorm;
99  int verbose;
101  double total_time;
102  double stage1_time;
103  double stage2_time;
104  double stage3_time;
105  double stage4_time;
107 #ifdef with_scalapack
108  int blacs_ctxt;
110 #endif
111 
112 } TridiagLU;
113 
114 int tridiagLU (double*,double*,double*,double*,int,int,void*,void*);
115 int tridiagLUGS (double*,double*,double*,double*,int,int,void*,void*);
116 int tridiagIterJacobi (double*,double*,double*,double*,int,int,void*,void*);
117 int tridiagLUInit (void*,void*);
118 
119 /* Block solvers */
120 int blocktridiagLU (double*,double*,double*,double*,int,int,int,void*,void*);
121 int blocktridiagIterJacobi (double*,double*,double*,double*,int,int,int,void*,void*);
122 
123 #ifdef with_scalapack
124 /* ScaLAPACK interface */
125 int tridiagScaLPK (double*,double*,double*,double*,int,int,void*,void*);
126 #endif
127 
128 #endif
double total_time
Definition: tridiagLU.h:101
int blocktridiagIterJacobi(double *, double *, double *, double *, int, int, int, void *, void *)
double stage4_time
Definition: tridiagLU.h:105
int tridiagLUInit(void *, void *)
Definition: tridiagLUInit.c:39
int verbose
Definition: tridiagLU.h:99
double stage3_time
Definition: tridiagLU.h:104
int tridiagLUGS(double *, double *, double *, double *, int, int, void *, void *)
Definition: tridiagLUGS.c:64
double exitnorm
Definition: tridiagLU.h:98
int tridiagIterJacobi(double *, double *, double *, double *, int, int, void *, void *)
int blocktridiagLU(double *, double *, double *, double *, int, int, int, void *, void *)
int evaluate_norm
Definition: tridiagLU.h:93
int exititer
Definition: tridiagLU.h:97
double stage2_time
Definition: tridiagLU.h:103
double stage1_time
Definition: tridiagLU.h:102
int tridiagLU(double *, double *, double *, double *, int, int, void *, void *)
Definition: tridiagLU.c:83
double rtol
Definition: tridiagLU.h:95
int maxiter
Definition: tridiagLU.h:94