xref: /petsc/src/ksp/pc/impls/tfs/tfs.h (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 #if !defined(__TFS_H)
3 #define __TFS_H
4 
5 /**********************************const.h*************************************
6 
7 Author: Henry M. Tufo III
8 
9 e-mail: hmt@cs.brown.edu
10 
11 snail-mail:
12 Division of Applied Mathematics
13 Brown University
14 Providence, RI 02912
15 
16 Last Modification:
17 6.21.97
18 ***********************************const.h************************************/
19 
20 /**********************************const.h*************************************
21 File Description:
22 -----------------
23 
24 ***********************************const.h************************************/
25 #include <petscsys.h>
26 #include <petscblaslapack.h>
27 
28 #define X  0
29 #define Y  1
30 #define Z  2
31 #define XY 3
32 #define XZ 4
33 #define YZ 5
34 
35 #define THRESH      0.2
36 #define N_HALF      4096
37 #define PRIV_BUF_SZ 45
38 
39 /*4096 8192 32768 65536 1048576 */
40 #define MAX_MSG_BUF 32768
41 
42 #define FULL    2
43 #define PARTIAL 1
44 #define NONE    0
45 
46 #define BYTE    8
47 #define BIT_0   0x1
48 #define BIT_1   0x2
49 #define BIT_2   0x4
50 #define BIT_3   0x8
51 #define BIT_4   0x10
52 #define BIT_5   0x20
53 #define BIT_6   0x40
54 #define BIT_7   0x80
55 #define TOP_BIT PETSC_MIN_INT
56 
57 #define C 0
58 
59 #define MAX_VEC     1674
60 #define FORMAT      30
61 #define MAX_COL_LEN 100
62 #define MAX_LINE    FORMAT *MAX_COL_LEN
63 #define DELIM       " \n \t"
64 #define LINE        12
65 #define C_LINE      80
66 
67 #define UT       5 /* dump upper 1/2 */
68 #define LT       6 /* dump lower 1/2 */
69 #define SYMM     8 /* we assume symm and dump upper 1/2 */
70 #define NON_SYMM 9
71 
72 #define ROW 10
73 #define COL 11
74 
75 #define EPS  1.0e-14
76 #define EPS2 1.0e-07
77 
78 #define MPI 1
79 #define NX  2
80 
81 #define LOG2(x) (PetscScalar) log((double)x) / log(2)
82 #define SWAP(a, b) \
83   temp = (a); \
84   (a)  = (b); \
85   (b)  = temp;
86 #define P_SWAP(a, b) \
87   ptr = (a); \
88   (a) = (b); \
89   (b) = ptr;
90 
91 #define MAX_FABS(x, y) (PetscAbsScalar(x) > PetscAbsScalar(y)) ? ((PetscScalar)x) : ((PetscScalar)y)
92 #define MIN_FABS(x, y) (PetscAbsScalar(x) < PetscAbsScalar(y)) ? ((PetscScalar)x) : ((PetscScalar)y)
93 
94 /* specer's existence ... can be done w/MAX_ABS */
95 #define EXISTS(x, y) ((x) == 0.0) ? (y) : (x)
96 
97 #define MULT_NEG_ONE(a) (a) *= -1;
98 #define NEG(a)          (a) |= BIT_31;
99 #define POS(a)          (a) &= INT_MAX;
100 
101 /**********************************types.h*************************************
102 
103 Author: Henry M. Tufo III
104 
105 e-mail: hmt@cs.brown.edu
106 
107 snail-mail:
108 Division of Applied Mathematics
109 Brown University
110 Providence, RI 02912
111 
112 Last Modification:
113 6.21.97
114 ***********************************types.h************************************/
115 
116 typedef PetscErrorCode (*vfp)(void *, void *, PetscInt, ...);
117 typedef PetscErrorCode (*rbfp)(PetscScalar *, PetscScalar *, PetscInt);
118 typedef PetscInt (*bfp)(void *, void *, PetscInt *, MPI_Datatype *);
119 
120 /***********************************comm.h*************************************
121 
122 Author: Henry M. Tufo III
123 
124 e-mail: hmt@cs.brown.edu
125 
126 snail-mail:
127 Division of Applied Mathematics
128 Brown University
129 Providence, RI 02912
130 
131 Last Modification:
132 6.21.97
133 ***********************************comm.h*************************************/
134 PETSC_INTERN PetscMPIInt PCTFS_my_id;
135 PETSC_INTERN PetscMPIInt PCTFS_num_nodes;
136 PETSC_INTERN PetscMPIInt PCTFS_floor_num_nodes;
137 PETSC_INTERN PetscMPIInt PCTFS_i_log2_num_nodes;
138 
139 PETSC_INTERN PetscErrorCode PCTFS_giop(PetscInt *, PetscInt *, PetscInt, PetscInt *);
140 PETSC_INTERN PetscErrorCode PCTFS_grop(PetscScalar *, PetscScalar *, PetscInt, PetscInt *);
141 PETSC_INTERN PetscErrorCode PCTFS_comm_init(void);
142 PETSC_INTERN PetscErrorCode PCTFS_giop_hc(PetscInt *, PetscInt *, PetscInt, PetscInt *, PetscInt);
143 PETSC_INTERN PetscErrorCode PCTFS_grop_hc(PetscScalar *, PetscScalar *, PetscInt, PetscInt *, PetscInt);
144 PETSC_INTERN PetscErrorCode PCTFS_ssgl_radd(PetscScalar *, PetscScalar *, PetscInt, PetscInt *);
145 
146 #define MSGTAG0 101
147 #define MSGTAG1 1001
148 #define MSGTAG2 76207
149 #define MSGTAG3 100001
150 #define MSGTAG4 163841
151 #define MSGTAG5 249439
152 #define MSGTAG6 10000001
153 
154 #define NON_UNIFORM 0
155 #define GL_MAX      1
156 #define GL_MIN      2
157 #define GL_MULT     3
158 #define GL_ADD      4
159 #define GL_B_XOR    5
160 #define GL_B_OR     6
161 #define GL_B_AND    7
162 #define GL_L_XOR    8
163 #define GL_L_OR     9
164 #define GL_L_AND    10
165 #define GL_MAX_ABS  11
166 #define GL_MIN_ABS  12
167 #define GL_EXISTS   13
168 
169 PETSC_INTERN PetscInt      *PCTFS_ivec_copy(PetscInt *, PetscInt *, PetscInt);
170 PETSC_INTERN PetscErrorCode PCTFS_ivec_zero(PetscInt *, PetscInt);
171 PETSC_INTERN PetscErrorCode PCTFS_ivec_set(PetscInt *, PetscInt, PetscInt);
172 
173 PETSC_INTERN PetscInt PCTFS_ivec_lb(PetscInt *, PetscInt);
174 PETSC_INTERN PetscInt PCTFS_ivec_ub(PetscInt *, PetscInt);
175 PETSC_INTERN PetscInt PCTFS_ivec_sum(PetscInt *, PetscInt);
176 PETSC_INTERN vfp      PCTFS_ivec_fct_addr(PetscInt);
177 
178 PETSC_INTERN PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *, PetscInt *, PetscInt, ...);
179 PETSC_INTERN PetscErrorCode PCTFS_ivec_max(PetscInt *, PetscInt *, PetscInt);
180 PETSC_INTERN PetscErrorCode PCTFS_ivec_min(PetscInt *, PetscInt *, PetscInt);
181 PETSC_INTERN PetscErrorCode PCTFS_ivec_mult(PetscInt *, PetscInt *, PetscInt);
182 PETSC_INTERN PetscErrorCode PCTFS_ivec_add(PetscInt *, PetscInt *, PetscInt);
183 PETSC_INTERN PetscErrorCode PCTFS_ivec_xor(PetscInt *, PetscInt *, PetscInt);
184 PETSC_INTERN PetscErrorCode PCTFS_ivec_or(PetscInt *, PetscInt *, PetscInt);
185 PETSC_INTERN PetscErrorCode PCTFS_ivec_and(PetscInt *, PetscInt *, PetscInt);
186 PETSC_INTERN PetscErrorCode PCTFS_ivec_lxor(PetscInt *, PetscInt *, PetscInt);
187 PETSC_INTERN PetscErrorCode PCTFS_ivec_lor(PetscInt *, PetscInt *, PetscInt);
188 PETSC_INTERN PetscErrorCode PCTFS_ivec_land(PetscInt *, PetscInt *, PetscInt);
189 PETSC_INTERN PetscErrorCode PCTFS_ivec_and3(PetscInt *, PetscInt *, PetscInt *, PetscInt);
190 
191 PETSC_INTERN PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *, PetscInt *, PetscInt);
192 PETSC_INTERN PetscErrorCode PCTFS_ivec_sort(PetscInt *, PetscInt);
193 PETSC_INTERN PetscErrorCode PCTFS_SMI_sort(void *, void *, PetscInt, PetscInt);
194 PETSC_INTERN PetscInt       PCTFS_ivec_binary_search(PetscInt, PetscInt *, PetscInt);
195 PETSC_INTERN PetscInt       PCTFS_ivec_linear_search(PetscInt, PetscInt *, PetscInt);
196 
197 PETSC_INTERN PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *, PetscInt **, PetscInt);
198 
199 #define SORT_INTEGER 1
200 #define SORT_INT_PTR 2
201 
202 PETSC_INTERN PetscErrorCode PCTFS_rvec_zero(PetscScalar *, PetscInt);
203 PETSC_INTERN PetscErrorCode PCTFS_rvec_one(PetscScalar *, PetscInt);
204 PETSC_INTERN PetscErrorCode PCTFS_rvec_set(PetscScalar *, PetscScalar, PetscInt);
205 PETSC_INTERN PetscErrorCode PCTFS_rvec_copy(PetscScalar *, PetscScalar *, PetscInt);
206 PETSC_INTERN PetscErrorCode PCTFS_rvec_scale(PetscScalar *, PetscScalar, PetscInt);
207 
208 PETSC_INTERN vfp            PCTFS_rvec_fct_addr(PetscInt);
209 PETSC_INTERN PetscErrorCode PCTFS_rvec_add(PetscScalar *, PetscScalar *, PetscInt);
210 PETSC_INTERN PetscErrorCode PCTFS_rvec_mult(PetscScalar *, PetscScalar *, PetscInt);
211 PETSC_INTERN PetscErrorCode PCTFS_rvec_max(PetscScalar *, PetscScalar *, PetscInt);
212 PETSC_INTERN PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *, PetscScalar *, PetscInt);
213 PETSC_INTERN PetscErrorCode PCTFS_rvec_min(PetscScalar *, PetscScalar *, PetscInt);
214 PETSC_INTERN PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *, PetscScalar *, PetscInt);
215 PETSC_INTERN PetscErrorCode PCTFS_vec_exists(PetscScalar *, PetscScalar *, PetscInt);
216 
217 /***********************************gs.h***************************************
218 
219 Author: Henry M. Tufo III
220 
221 e-mail: hmt@cs.brown.edu
222 
223 snail-mail:
224 Division of Applied Mathematics
225 Brown University
226 Providence, RI 02912
227 
228 Last Modification:
229 6.21.97
230 ************************************gs.h**************************************/
231 
232 typedef struct gather_scatter_id *PCTFS_gs_ADT;
233 
234 PETSC_INTERN PCTFS_gs_ADT   PCTFS_gs_init(PetscInt *, PetscInt, PetscInt);
235 PETSC_INTERN PetscErrorCode PCTFS_gs_gop_vec(PCTFS_gs_ADT, PetscScalar *, const char *, PetscInt);
236 PETSC_INTERN PetscErrorCode PCTFS_gs_gop_hc(PCTFS_gs_ADT, PetscScalar *, const char *, PetscInt);
237 PETSC_INTERN PetscErrorCode PCTFS_gs_free(PCTFS_gs_ADT);
238 PETSC_INTERN PetscErrorCode PCTFS_gs_init_msg_buf_sz(PetscInt);
239 PETSC_INTERN PetscErrorCode PCTFS_gs_init_vec_sz(PetscInt);
240 
241 /*************************************xxt.h************************************
242 Module Name: xxt
243 Module Info: need xxt.{c,h} gs.{c,h} comm.{c,h} ivec.{c,h} error.{c,h}
244 
245 author:  Henry M. Tufo III
246 e-mail:  hmt@asci.uchicago.edu
247 contact:
248 +--------------------------------+--------------------------------+
249 |MCS Division - Building 221     |Department of Computer Science  |
250 |Argonne National Laboratory     |Ryerson 152                     |
251 |9700 S. Cass Avenue             |The University of Chicago       |
252 |Argonne, IL  60439              |Chicago, IL  60637              |
253 |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
254 +--------------------------------+--------------------------------+
255 
256 Last Modification: 3.20.01
257 **************************************xxt.h***********************************/
258 
259 typedef struct xxt_CDT *xxt_ADT;
260 
261 /*************************************xxt.h************************************
262 Function: XXT_new()
263 
264 Return: ADT ptr or NULL upon failure.
265 Description: This function allocates and returns an xxt handle
266 Usage: xxt_handle = xxt_new();
267 **************************************xxt.h***********************************/
268 PETSC_INTERN xxt_ADT XXT_new(void);
269 
270 /*************************************xxt.h************************************
271 Function: XXT_free()
272 
273 Input : pointer to ADT.
274 
275 Description: This function frees the storage associated with an xxt handle
276 Usage: XXT_free(xxt_handle);
277 **************************************xxt.h***********************************/
278 PETSC_INTERN PetscInt XXT_free(xxt_ADT);
279 
280 /*************************************xxt.h************************************
281 Function: XXT_factor
282 
283 Input : ADT ptr,  and pointer to object
284 Return: 0 on failure, 1 on success
285 Description: This function sets the xxt solver
286 
287 xxt assumptions: given n rows of global coarse matrix (E_loc) where
288    o global dofs N = sum_p(n), p=0,P-1
289    (i.e. row dist. with no dof replication)
290    (5.21.00 will handle dif replication case)
291    o m is the number of columns in E_loc (m>=n)
292    o local2global holds global number of column i (i=0,...,m-1)
293    o local2global holds global number of row    i (i=0,...,n-1)
294    o mylocmatvec performs E_loc . x_loc where x_loc is an vector of
295    length m in 1-1 correspondence with local2global
296    (note that gs package takes care of communication).
297    (note do not zero out upper m-n entries!)
298    o mylocmatvec(void *grid_data, double *in, double *out)
299 
300 ML beliefs/usage: move this to to ML_XXT_factor routine
301    o my_ml holds address of ML struct associated w/E_loc, grid_data, grid_tag
302    o grid_tag, grid_data, my_ml used in
303       ML_Set_CSolve(my_ml, grid_tag, grid_data, ML_Do_CoarseDirect);
304    o grid_data used in
305       A_matvec(grid_data,v,u);
306 
307 Usage:
308 **************************************xxt.h***********************************/
309 PETSC_INTERN PetscErrorCode XXT_factor(xxt_ADT,                                                  /* prev. allocated xxt  handle */
310                                        PetscInt *,                                               /* global column mapping       */
311                                        PetscInt,                                                 /* local num rows              */
312                                        PetscInt,                                                 /* local num cols              */
313                                        PetscErrorCode (*)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc         */
314                                        void *);                                                  /* grid data for matvec        */
315 
316 /*************************************xxt.h************************************
317 Function: XXT_solve
318 
319 Input : ADT ptr, b (rhs)
320 Output: x (soln)
321 Return:
322 Description: This function performs x = E^-1.b
323 Usage:
324 XXT_solve(xxt_handle, double *x, double *b)
325 XXT_solve(xxt_handle, double *x, NULL)
326 assumes x has been initialized to be b
327 **************************************xxt.h***********************************/
328 PETSC_INTERN PetscErrorCode XXT_solve(xxt_ADT, PetscScalar *, PetscScalar *);
329 
330 /*************************************xxt.h************************************
331 Function: XXT_stats
332 
333 Input : handle
334 **************************************xxt.h***********************************/
335 PETSC_INTERN PetscErrorCode XXT_stats(xxt_ADT);
336 
337 /*************************************xxt.h************************************
338 Function: XXT_sp_1()
339 
340 Input : pointer to ADT
341 Output:
342 Return:
343 Description: sets xxt parameter 1 in xxt_handle
344 Usage: implement later
345 
346 void XXT_sp_1(xxt_handle,parameter 1 value)
347 **************************************xxt.h***********************************/
348 
349 /*************************************xyt.h************************************
350 Module Name: xyt
351 Module Info: need xyt.{c,h} gs.{c,h} comm.{c,h} ivec.{c,h} error.{c,h}
352 
353 author:  Henry M. Tufo III
354 e-mail:  hmt@asci.uchicago.edu
355 contact:
356 +--------------------------------+--------------------------------+
357 |MCS Division - Building 221     |Department of Computer Science  |
358 |Argonne National Laboratory     |Ryerson 152                     |
359 |9700 S. Cass Avenue             |The University of Chicago       |
360 |Argonne, IL  60439              |Chicago, IL  60637              |
361 |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
362 +--------------------------------+--------------------------------+
363 
364 Last Modification: 3.20.01
365 **************************************xyt.h***********************************/
366 
367 typedef struct xyt_CDT *xyt_ADT;
368 
369 /*************************************xyt.h************************************
370 Function: XYT_new()
371 
372 Return: ADT ptr or NULL upon failure.
373 Description: This function allocates and returns an xyt handle
374 Usage: xyt_handle = xyt_new();
375 **************************************xyt.h***********************************/
376 PETSC_INTERN xyt_ADT XYT_new(void);
377 
378 /*************************************xyt.h************************************
379 Function: XYT_free()
380 
381 Input : pointer to ADT.
382 Description: This function frees the storage associated with an xyt handle
383 Usage: XYT_free(xyt_handle);
384 **************************************xyt.h***********************************/
385 PETSC_INTERN PetscErrorCode XYT_free(xyt_ADT);
386 
387 /*************************************xyt.h************************************
388 Function: XYT_factor
389 
390 Input : ADT ptr,  and pointer to object
391 Output:
392 Return: 0 on failure, 1 on success
393 Description: This function sets the xyt solver
394 
395 xyt assumptions: given n rows of global coarse matrix (E_loc) where
396    o global dofs N = sum_p(n), p=0,P-1
397    (i.e. row dist. with no dof replication)
398    (5.21.00 will handle dif replication case)
399    o m is the number of columns in E_loc (m>=n)
400    o local2global holds global number of column i (i=0,...,m-1)
401    o local2global holds global number of row    i (i=0,...,n-1)
402    o mylocmatvec performs E_loc . x_loc where x_loc is an vector of
403    length m in 1-1 correspondence with local2global
404    (note that gs package takes care of communication).
405    (note do not zero out upper m-n entries!)
406    o mylocmatvec(void *grid_data, double *in, double *out)
407 
408 ML beliefs/usage: move this to to ML_XYT_factor routine
409    o my_ml holds address of ML struct associated w/E_loc, grid_data, grid_tag
410    o grid_tag, grid_data, my_ml used in
411       ML_Set_CSolve(my_ml, grid_tag, grid_data, ML_Do_CoarseDirect);
412    o grid_data used in
413       A_matvec(grid_data,v,u);
414 
415 Usage:
416 **************************************xyt.h***********************************/
417 PETSC_INTERN PetscErrorCode XYT_factor(xyt_ADT,                                                  /* prev. allocated xyt  handle */
418                                        PetscInt *,                                               /* global column mapping       */
419                                        PetscInt,                                                 /* local num rows              */
420                                        PetscInt,                                                 /* local num cols              */
421                                        PetscErrorCode (*)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc         */
422                                        void *);                                                  /* grid data for matvec        */
423 
424 /*************************************xyt.h************************************
425 Function: XYT_solve
426 
427 Input : ADT ptr, b (rhs)
428 Output: x (soln)
429 Return:
430 Description: This function performs x = E^-1.b
431 Usage: XYT_solve(xyt_handle, double *x, double *b)
432 **************************************xyt.h***********************************/
433 PETSC_INTERN PetscErrorCode XYT_solve(xyt_ADT, PetscScalar *, PetscScalar *);
434 
435 /*************************************xyt.h************************************
436 Function: XYT_stats
437 
438 Input : handle
439 **************************************xyt.h***********************************/
440 PETSC_INTERN PetscErrorCode XYT_stats(xyt_ADT);
441 
442 /********************************bit_mask.h************************************
443 
444 Author: Henry M. Tufo III
445 
446 e-mail: hmt@cs.brown.edu
447 
448 snail-mail:
449 Division of Applied Mathematics
450 Brown University
451 Providence, RI 02912
452 
453 Last Modification:
454 11.21.97
455 *********************************bit_mask.h***********************************/
456 PETSC_INTERN PetscInt       PCTFS_div_ceil(PetscInt, PetscInt);
457 PETSC_INTERN PetscErrorCode PCTFS_set_bit_mask(PetscInt *, PetscInt, PetscInt);
458 PETSC_INTERN PetscInt       PCTFS_len_bit_mask(PetscInt);
459 PETSC_INTERN PetscInt       PCTFS_ct_bits(char *, PetscInt);
460 PETSC_INTERN PetscErrorCode PCTFS_bm_to_proc(char *, PetscInt, PetscInt *);
461 PETSC_INTERN PetscInt       PCTFS_len_buf(PetscInt, PetscInt);
462 
463 #endif
464