Actual source code: tfs.h
petsc-3.11.4 2019-09-28
5: /**********************************const.h*************************************
7: Author: Henry M. Tufo III
9: e-mail: hmt@cs.brown.edu
11: snail-mail:
12: Division of Applied Mathematics
13: Brown University
14: Providence, RI 02912
16: Last Modification:
17: 6.21.97
18: ***********************************const.h************************************/
20: /**********************************const.h*************************************
21: File Description:
22: -----------------
24: ***********************************const.h************************************/
25: #include <petscsys.h>
26: #include <petscblaslapack.h>
28: #define X 0
29: #define Y 1
30: #define Z 2
31: #define XY 3
32: #define XZ 4
33: #define YZ 5
36: #define THRESH 0.2
37: #define N_HALF 4096
38: #define PRIV_BUF_SZ 45
40: /*4096 8192 32768 65536 1048576 */
41: #define MAX_MSG_BUF 32768
43: #define FULL 2
44: #define PARTIAL 1
45: #define NONE 0
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
58: #define C 0
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
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
74: #define ROW 10
75: #define COL 11
77: #define EPS 1.0e-14
78: #define EPS2 1.0e-07
81: #define MPI 1
82: #define NX 2
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;
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)
91: /* specer's existence ... can be done w/MAX_ABS */
92: #define EXISTS(x,y) ((x)==0.0) ? (y) : (x)
94: #define MULT_NEG_ONE(a) (a) *= -1;
95: #define NEG(a) (a) |= BIT_31;
96: #define POS(a) (a) &= INT_MAX;
101: /**********************************types.h*************************************
103: Author: Henry M. Tufo III
105: e-mail: hmt@cs.brown.edu
107: snail-mail:
108: Division of Applied Mathematics
109: Brown University
110: Providence, RI 02912
112: Last Modification:
113: 6.21.97
114: ***********************************types.h************************************/
116: typedef PetscErrorCode (*vfp)(void*,void*,PetscInt,...);
117: typedef PetscErrorCode (*rbfp)(PetscScalar*, PetscScalar*, PetscInt);
118: typedef PetscInt (*bfp)(void*, void*, PetscInt*, MPI_Datatype*);
120: /***********************************comm.h*************************************
122: Author: Henry M. Tufo III
124: e-mail: hmt@cs.brown.edu
126: snail-mail:
127: Division of Applied Mathematics
128: Brown University
129: Providence, RI 02912
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;
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*);
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
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
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);
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);
178: PETSC_INTERN PetscErrorCode PCTFS_ivec_non_uniform(PetscInt*,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);
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);
197: PETSC_INTERN PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt*,PetscInt**,PetscInt);
199: #define SORT_INTEGER 1
200: #define SORT_INT_PTR 2
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);
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);
217: /***********************************gs.h***************************************
219: Author: Henry M. Tufo III
221: e-mail: hmt@cs.brown.edu
223: snail-mail:
224: Division of Applied Mathematics
225: Brown University
226: Providence, RI 02912
228: Last Modification:
229: 6.21.97
230: ************************************gs.h**************************************/
232: typedef struct gather_scatter_id *PCTFS_gs_ADT;
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);
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}
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: +--------------------------------+--------------------------------+
256: Last Modification: 3.20.01
257: **************************************xxt.h***********************************/
259: typedef struct xxt_CDT *xxt_ADT;
262: /*************************************xxt.h************************************
263: Function: XXT_new()
265: Return: ADT ptr or NULL upon failure.
266: Description: This function allocates and returns an xxt handle
267: Usage: xxt_handle = xxt_new();
268: **************************************xxt.h***********************************/
269: PETSC_INTERN xxt_ADT XXT_new(void);
272: /*************************************xxt.h************************************
273: Function: XXT_free()
275: Input : pointer to ADT.
277: Description: This function frees the storage associated with an xxt handle
278: Usage: XXT_free(xxt_handle);
279: **************************************xxt.h***********************************/
280: PETSC_INTERN PetscInt XXT_free(xxt_ADT);
283: /*************************************xxt.h************************************
284: Function: XXT_factor
286: Input : ADT ptr, and pointer to object
287: Return: 0 on failure, 1 on success
288: Description: This function sets the xxt solver
290: xxt assumptions: given n rows of global coarse matrix (E_loc) where
291: o global dofs N = sum_p(n), p=0,P-1
292: (i.e. row dist. with no dof replication)
293: (5.21.00 will handle dif replication case)
294: o m is the number of columns in E_loc (m>=n)
295: o local2global holds global number of column i (i=0,...,m-1)
296: o local2global holds global number of row i (i=0,...,n-1)
297: o mylocmatvec performs E_loc . x_loc where x_loc is an vector of
298: length m in 1-1 correspondence with local2global
299: (note that gs package takes care of communication).
300: (note do not zero out upper m-n entries!)
301: o mylocmatvec(void *grid_data, double *in, double *out)
303: ML beliefs/usage: move this to to ML_XXT_factor routine
304: o my_ml holds address of ML struct associated w/E_loc, grid_data, grid_tag
305: o grid_tag, grid_data, my_ml used in
306: ML_Set_CSolve(my_ml, grid_tag, grid_data, ML_Do_CoarseDirect);
307: o grid_data used in
308: A_matvec(grid_data,v,u);
310: Usage:
311: **************************************xxt.h***********************************/
312: PETSC_INTERN PetscInt XXT_factor(xxt_ADT, /* prev. allocated xxt handle */
313: PetscInt*, /* global column mapping */
314: PetscInt, /* local num rows */
315: PetscInt, /* local num cols */
316: PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */
317: void*); /* grid data for matvec */
320: /*************************************xxt.h************************************
321: Function: XXT_solve
323: Input : ADT ptr, b (rhs)
324: Output: x (soln)
325: Return:
326: Description: This function performs x = E^-1.b
327: Usage:
328: XXT_solve(xxt_handle, double *x, double *b)
329: XXT_solve(xxt_handle, double *x, NULL)
330: assumes x has been initialized to be b
331: **************************************xxt.h***********************************/
332: PETSC_INTERN PetscInt XXT_solve(xxt_ADT,PetscScalar*,PetscScalar*);
334: /*************************************xxt.h************************************
335: Function: XXT_stats
337: Input : handle
338: **************************************xxt.h***********************************/
339: PETSC_INTERN PetscInt XXT_stats(xxt_ADT);
342: /*************************************xxt.h************************************
343: Function: XXT_sp_1()
345: Input : pointer to ADT
346: Output:
347: Return:
348: Description: sets xxt parameter 1 in xxt_handle
349: Usage: implement later
351: void XXT_sp_1(xxt_handle,parameter 1 value)
352: **************************************xxt.h***********************************/
355: /*************************************xyt.h************************************
356: Module Name: xyt
357: Module Info: need xyt.{c,h} gs.{c,h} comm.{c,h} ivec.{c,h} error.{c,h}
359: author: Henry M. Tufo III
360: e-mail: hmt@asci.uchicago.edu
361: contact:
362: +--------------------------------+--------------------------------+
363: |MCS Division - Building 221 |Department of Computer Science |
364: |Argonne National Laboratory |Ryerson 152 |
365: |9700 S. Cass Avenue |The University of Chicago |
366: |Argonne, IL 60439 |Chicago, IL 60637 |
367: |(630) 252-5354/5986 ph/fx |(773) 702-6019/8487 ph/fx |
368: +--------------------------------+--------------------------------+
370: Last Modification: 3.20.01
371: **************************************xyt.h***********************************/
373: typedef struct xyt_CDT *xyt_ADT;
376: /*************************************xyt.h************************************
377: Function: XYT_new()
379: Return: ADT ptr or NULL upon failure.
380: Description: This function allocates and returns an xyt handle
381: Usage: xyt_handle = xyt_new();
382: **************************************xyt.h***********************************/
383: PETSC_INTERN xyt_ADT XYT_new(void);
386: /*************************************xyt.h************************************
387: Function: XYT_free()
389: Input : pointer to ADT.
390: Description: This function frees the storage associated with an xyt handle
391: Usage: XYT_free(xyt_handle);
392: **************************************xyt.h***********************************/
393: PETSC_INTERN PetscInt XYT_free(xyt_ADT);
396: /*************************************xyt.h************************************
397: Function: XYT_factor
399: Input : ADT ptr, and pointer to object
400: Output:
401: Return: 0 on failure, 1 on success
402: Description: This function sets the xyt solver
404: xyt assumptions: given n rows of global coarse matrix (E_loc) where
405: o global dofs N = sum_p(n), p=0,P-1
406: (i.e. row dist. with no dof replication)
407: (5.21.00 will handle dif replication case)
408: o m is the number of columns in E_loc (m>=n)
409: o local2global holds global number of column i (i=0,...,m-1)
410: o local2global holds global number of row i (i=0,...,n-1)
411: o mylocmatvec performs E_loc . x_loc where x_loc is an vector of
412: length m in 1-1 correspondence with local2global
413: (note that gs package takes care of communication).
414: (note do not zero out upper m-n entries!)
415: o mylocmatvec(void *grid_data, double *in, double *out)
417: ML beliefs/usage: move this to to ML_XYT_factor routine
418: o my_ml holds address of ML struct associated w/E_loc, grid_data, grid_tag
419: o grid_tag, grid_data, my_ml used in
420: ML_Set_CSolve(my_ml, grid_tag, grid_data, ML_Do_CoarseDirect);
421: o grid_data used in
422: A_matvec(grid_data,v,u);
424: Usage:
425: **************************************xyt.h***********************************/
426: PETSC_INTERN PetscInt XYT_factor(xyt_ADT, /* prev. allocated xyt handle */
427: PetscInt*, /* global column mapping */
428: PetscInt, /* local num rows */
429: PetscInt, /* local num cols */
430: PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */
431: void*); /* grid data for matvec */
434: /*************************************xyt.h************************************
435: Function: XYT_solve
437: Input : ADT ptr, b (rhs)
438: Output: x (soln)
439: Return:
440: Description: This function performs x = E^-1.b
441: Usage: XYT_solve(xyt_handle, double *x, double *b)
442: **************************************xyt.h***********************************/
443: PETSC_INTERN PetscInt XYT_solve(xyt_ADT,PetscScalar*,PetscScalar*);
446: /*************************************xyt.h************************************
447: Function: XYT_stats
449: Input : handle
450: **************************************xyt.h***********************************/
451: PETSC_INTERN PetscInt XYT_stats(xyt_ADT);
454: /********************************bit_mask.h************************************
456: Author: Henry M. Tufo III
458: e-mail: hmt@cs.brown.edu
460: snail-mail:
461: Division of Applied Mathematics
462: Brown University
463: Providence, RI 02912
465: Last Modification:
466: 11.21.97
467: *********************************bit_mask.h***********************************/
468: PETSC_INTERN PetscInt PCTFS_div_ceil(PetscInt,PetscInt);
469: PETSC_INTERN PetscErrorCode PCTFS_set_bit_mask(PetscInt*,PetscInt,PetscInt);
470: PETSC_INTERN PetscInt PCTFS_len_bit_mask(PetscInt);
471: PETSC_INTERN PetscInt PCTFS_ct_bits(char*, PetscInt);
472: PETSC_INTERN PetscErrorCode PCTFS_bm_to_proc(char*,PetscInt,PetscInt*);
473: PETSC_INTERN PetscInt PCTFS_len_buf(PetscInt, PetscInt);
475: #endif