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