Actual source code: tfs.h

  1: #pragma once
  2: #if defined(__GNUC__) || defined(__clang__)
  3:   #pragma GCC diagnostic ignored "-Wconversion"
  4: #endif

  6: /**********************************const.h*************************************

  8: Author: Henry M. Tufo III

 10: e-mail: hmt@cs.brown.edu

 12: snail-mail:
 13: Division of Applied Mathematics
 14: Brown University
 15: Providence, RI 02912

 17: Last Modification:
 18: 6.21.97
 19: ***********************************const.h************************************/

 21: /**********************************const.h*************************************
 22: File Description:
 23: -----------------

 25: ***********************************const.h************************************/
 26: #include <petscsys.h>
 27: #include <petscblaslapack.h>

 29: #define X  0
 30: #define Y  1
 31: #define Z  2
 32: #define XY 3
 33: #define XZ 4
 34: #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_INT_MIN

 58: #define C 0

 60: #define MAX_VEC     1674
 61: #define FORMAT      30
 62: #define MAX_COL_LEN 100
 63: #define MAX_LINE    FORMAT *MAX_COL_LEN
 64: #define DELIM       " \n \t"
 65: #define LINE        12
 66: #define C_LINE      80

 68: #define UT       5 /* dump upper 1/2 */
 69: #define LT       6 /* dump lower 1/2 */
 70: #define SYMM     8 /* we assume symm and dump upper 1/2 */
 71: #define NON_SYMM 9

 73: #define ROW 10
 74: #define COL 11

 76: #define EPS  1.0e-14
 77: #define EPS2 1.0e-07

 79: #define MPI 1
 80: #define NX  2

 82: #define LOG2(x) (PetscScalar)log((double)x) / log(2)
 83: #define SWAP(a, b) \
 84:   temp = (a); \
 85:   (a)  = (b); \
 86:   (b)  = temp;
 87: #define P_SWAP(a, b) \
 88:   ptr = (a); \
 89:   (a) = (b); \
 90:   (b) = ptr;

 92: #define MAX_FABS(x, y) (PetscAbsScalar(x) > PetscAbsScalar(y)) ? ((PetscScalar)x) : ((PetscScalar)y)
 93: #define MIN_FABS(x, y) (PetscAbsScalar(x) < PetscAbsScalar(y)) ? ((PetscScalar)x) : ((PetscScalar)y)

 95: /* specer's existence ... can be done w/MAX_ABS */
 96: #define EXISTS(x, y) ((x) == 0.0) ? (y) : (x)

 98: #define MULT_NEG_ONE(a) (a) *= -1;
 99: #define NEG(a)          (a) |= BIT_31;
100: #define POS(a)          (a) &= INT_MAX;

102: /**********************************types.h*************************************

104: Author: Henry M. Tufo III

106: e-mail: hmt@cs.brown.edu

108: snail-mail:
109: Division of Applied Mathematics
110: Brown University
111: Providence, RI 02912

113: Last Modification:
114: 6.21.97
115: ***********************************types.h************************************/

117: typedef PetscErrorCode (*vfp)(void *, void *, PetscInt, ...);
118: typedef PetscErrorCode (*rbfp)(PetscScalar *, PetscScalar *, PetscInt);
119: typedef PetscInt (*bfp)(void *, void *, PetscInt *, MPI_Datatype *);

121: /***********************************comm.h*************************************

123: Author: Henry M. Tufo III

125: e-mail: hmt@cs.brown.edu

127: snail-mail:
128: Division of Applied Mathematics
129: Brown University
130: Providence, RI 02912

132: Last Modification:
133: 6.21.97
134: ***********************************comm.h*************************************/
135: PETSC_INTERN PetscMPIInt PCTFS_my_id;
136: PETSC_INTERN PetscMPIInt PCTFS_num_nodes;
137: PETSC_INTERN PetscMPIInt PCTFS_floor_num_nodes;
138: PETSC_INTERN PetscMPIInt PCTFS_i_log2_num_nodes;

140: PETSC_INTERN PetscErrorCode PCTFS_giop(PetscInt *, PetscInt *, PetscInt, PetscInt *);
141: PETSC_INTERN PetscErrorCode PCTFS_grop(PetscScalar *, PetscScalar *, PetscInt, PetscInt *);
142: PETSC_INTERN PetscErrorCode PCTFS_comm_init(void);
143: PETSC_INTERN PetscErrorCode PCTFS_giop_hc(PetscInt *, PetscInt *, PetscInt, PetscInt *, PetscInt);
144: PETSC_INTERN PetscErrorCode PCTFS_grop_hc(PetscScalar *, PetscScalar *, PetscInt, PetscInt *, PetscInt);
145: PETSC_INTERN PetscErrorCode PCTFS_ssgl_radd(PetscScalar *, PetscScalar *, PetscInt, PetscInt *);

147: #define MSGTAG0 101
148: #define MSGTAG1 1001
149: #define MSGTAG2 76207
150: #define MSGTAG3 100001
151: #define MSGTAG4 163841
152: #define MSGTAG5 249439
153: #define MSGTAG6 10000001

155: #define NON_UNIFORM 0
156: #define GL_MAX      1
157: #define GL_MIN      2
158: #define GL_MULT     3
159: #define GL_ADD      4
160: #define GL_B_XOR    5
161: #define GL_B_OR     6
162: #define GL_B_AND    7
163: #define GL_L_XOR    8
164: #define GL_L_OR     9
165: #define GL_L_AND    10
166: #define GL_MAX_ABS  11
167: #define GL_MIN_ABS  12
168: #define GL_EXISTS   13

170: PETSC_INTERN PetscInt      *PCTFS_ivec_copy(PetscInt *, PetscInt *, PetscInt);
171: PETSC_INTERN PetscErrorCode PCTFS_ivec_zero(PetscInt *, PetscInt);
172: PETSC_INTERN PetscErrorCode PCTFS_ivec_set(PetscInt *, PetscInt, PetscInt);

174: PETSC_INTERN PetscInt PCTFS_ivec_lb(PetscInt *, PetscInt);
175: PETSC_INTERN PetscInt PCTFS_ivec_ub(PetscInt *, PetscInt);
176: PETSC_INTERN PetscInt PCTFS_ivec_sum(PetscInt *, PetscInt);
177: PETSC_INTERN vfp      PCTFS_ivec_fct_addr(PetscInt);

179: PETSC_INTERN PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *, PetscInt *, PetscInt, ...);
180: PETSC_INTERN PetscErrorCode PCTFS_ivec_max(PetscInt *, PetscInt *, PetscInt);
181: PETSC_INTERN PetscErrorCode PCTFS_ivec_min(PetscInt *, PetscInt *, PetscInt);
182: PETSC_INTERN PetscErrorCode PCTFS_ivec_mult(PetscInt *, PetscInt *, PetscInt);
183: PETSC_INTERN PetscErrorCode PCTFS_ivec_add(PetscInt *, PetscInt *, PetscInt);
184: PETSC_INTERN PetscErrorCode PCTFS_ivec_xor(PetscInt *, PetscInt *, PetscInt);
185: PETSC_INTERN PetscErrorCode PCTFS_ivec_or(PetscInt *, PetscInt *, PetscInt);
186: PETSC_INTERN PetscErrorCode PCTFS_ivec_and(PetscInt *, PetscInt *, PetscInt);
187: PETSC_INTERN PetscErrorCode PCTFS_ivec_lxor(PetscInt *, PetscInt *, PetscInt);
188: PETSC_INTERN PetscErrorCode PCTFS_ivec_lor(PetscInt *, PetscInt *, PetscInt);
189: PETSC_INTERN PetscErrorCode PCTFS_ivec_land(PetscInt *, PetscInt *, PetscInt);
190: PETSC_INTERN PetscErrorCode PCTFS_ivec_and3(PetscInt *, PetscInt *, PetscInt *, PetscInt);

192: PETSC_INTERN PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *, PetscInt *, PetscInt);
193: PETSC_INTERN PetscErrorCode PCTFS_ivec_sort(PetscInt *, PetscInt);
194: PETSC_INTERN PetscErrorCode PCTFS_SMI_sort(void *, void *, PetscInt, PetscInt);
195: PETSC_INTERN PetscInt       PCTFS_ivec_binary_search(PetscInt, PetscInt *, PetscInt);
196: PETSC_INTERN PetscInt       PCTFS_ivec_linear_search(PetscInt, PetscInt *, PetscInt);

198: PETSC_INTERN PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *, PetscInt **, PetscInt);

200: #define SORT_INTEGER 1
201: #define SORT_INT_PTR 2

203: PETSC_INTERN PetscErrorCode PCTFS_rvec_zero(PetscScalar *, PetscInt);
204: PETSC_INTERN PetscErrorCode PCTFS_rvec_one(PetscScalar *, PetscInt);
205: PETSC_INTERN PetscErrorCode PCTFS_rvec_set(PetscScalar *, PetscScalar, PetscInt);
206: PETSC_INTERN PetscErrorCode PCTFS_rvec_copy(PetscScalar *, PetscScalar *, PetscInt);
207: PETSC_INTERN PetscErrorCode PCTFS_rvec_scale(PetscScalar *, PetscScalar, PetscInt);

209: PETSC_INTERN vfp            PCTFS_rvec_fct_addr(PetscInt);
210: PETSC_INTERN PetscErrorCode PCTFS_rvec_add(PetscScalar *, PetscScalar *, PetscInt);
211: PETSC_INTERN PetscErrorCode PCTFS_rvec_mult(PetscScalar *, PetscScalar *, PetscInt);
212: PETSC_INTERN PetscErrorCode PCTFS_rvec_max(PetscScalar *, PetscScalar *, PetscInt);
213: PETSC_INTERN PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *, PetscScalar *, PetscInt);
214: PETSC_INTERN PetscErrorCode PCTFS_rvec_min(PetscScalar *, PetscScalar *, PetscInt);
215: PETSC_INTERN PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *, PetscScalar *, PetscInt);
216: PETSC_INTERN PetscErrorCode PCTFS_vec_exists(PetscScalar *, PetscScalar *, PetscInt);

218: /***********************************gs.h***************************************

220: Author: Henry M. Tufo III

222: e-mail: hmt@cs.brown.edu

224: snail-mail:
225: Division of Applied Mathematics
226: Brown University
227: Providence, RI 02912

229: Last Modification:
230: 6.21.97
231: ************************************gs.h**************************************/

233: typedef struct gather_scatter_id *PCTFS_gs_ADT;

235: PETSC_INTERN PCTFS_gs_ADT   PCTFS_gs_init(PetscInt *, PetscInt, PetscInt);
236: PETSC_INTERN PetscErrorCode PCTFS_gs_gop_vec(PCTFS_gs_ADT, PetscScalar *, const char *, PetscInt);
237: PETSC_INTERN PetscErrorCode PCTFS_gs_gop_hc(PCTFS_gs_ADT, PetscScalar *, const char *, PetscInt);
238: PETSC_INTERN PetscErrorCode PCTFS_gs_free(PCTFS_gs_ADT);
239: PETSC_INTERN PetscErrorCode PCTFS_gs_init_msg_buf_sz(PetscInt);
240: PETSC_INTERN PetscErrorCode PCTFS_gs_init_vec_sz(PetscInt);

242: /*************************************xxt.h************************************
243: Module Name: xxt
244: Module Info: need xxt.{c,h} gs.{c,h} comm.{c,h} ivec.{c,h} error.{c,h}

246: author:  Henry M. Tufo III
247: e-mail:  hmt@asci.uchicago.edu
248: contact:
249: +--------------------------------+--------------------------------+
250: |MCS Division - Building 221     |Department of Computer Science  |
251: |Argonne National Laboratory     |Ryerson 152                     |
252: |9700 S. Cass Avenue             |The University of Chicago       |
253: |Argonne, IL  60439              |Chicago, IL  60637              |
254: |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
255: +--------------------------------+--------------------------------+

257: Last Modification: 3.20.01
258: **************************************xxt.h***********************************/

260: 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);

271: /*************************************xxt.h************************************
272: Function: XXT_free()

274: Input : pointer to ADT.

276: Description: This function frees the storage associated with an xxt handle
277: Usage: XXT_free(xxt_handle);
278: **************************************xxt.h***********************************/
279: PETSC_INTERN PetscErrorCode XXT_free(xxt_ADT);

281: /*************************************xxt.h************************************
282: Function: XXT_factor

284: Input : ADT ptr,  and pointer to object
285: Return: 0 on failure, 1 on success
286: Description: This function sets the xxt solver

288: xxt assumptions: given n rows of global coarse matrix (E_loc) where
289:    o global dofs N = sum_p(n), p=0,P-1
290:    (i.e. row dist. with no dof replication)
291:    (5.21.00 will handle dif replication case)
292:    o m is the number of columns in E_loc (m>=n)
293:    o local2global holds global number of column i (i=0,...,m-1)
294:    o local2global holds global number of row    i (i=0,...,n-1)
295:    o mylocmatvec performs E_loc . x_loc where x_loc is an vector of
296:    length m in 1-1 correspondence with local2global
297:    (note that gs package takes care of communication).
298:    (note do not zero out upper m-n entries!)
299:    o mylocmatvec(void *grid_data, double *in, double *out)

301: ML beliefs/usage: move this to ML_XXT_factor routine
302:    o my_ml holds address of ML struct associated w/E_loc, grid_data, grid_tag
303:    o grid_tag, grid_data, my_ml used in
304:       ML_Set_CSolve(my_ml, grid_tag, grid_data, ML_Do_CoarseDirect);
305:    o grid_data used in
306:       A_matvec(grid_data,v,u);

308: Usage:
309: **************************************xxt.h***********************************/
310: PETSC_INTERN PetscErrorCode XXT_factor(xxt_ADT,                                                  /* prev. allocated xxt  handle */
311:                                        PetscInt *,                                               /* global column mapping       */
312:                                        PetscInt,                                                 /* local num rows              */
313:                                        PetscInt,                                                 /* local num cols              */
314:                                        PetscErrorCode (*)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc         */
315:                                        void *);                                                  /* grid data for matvec        */

317: /*************************************xxt.h************************************
318: Function: XXT_solve

320: Input : ADT ptr, b (rhs)
321: Output: x (soln)
322: Return:
323: Description: This function performs x = E^-1.b
324: Usage:
325: XXT_solve(xxt_handle, double *x, double *b)
326: XXT_solve(xxt_handle, double *x, NULL)
327: assumes x has been initialized to be b
328: **************************************xxt.h***********************************/
329: PETSC_INTERN PetscErrorCode XXT_solve(xxt_ADT, PetscScalar *, PetscScalar *);

331: /*************************************xxt.h************************************
332: Function: XXT_sp_1()

334: Input : pointer to ADT
335: Output:
336: Return:
337: Description: sets xxt parameter 1 in xxt_handle
338: Usage: implement later

340: void XXT_sp_1(xxt_handle,parameter 1 value)
341: **************************************xxt.h***********************************/

343: /*************************************xyt.h************************************
344: Module Name: xyt
345: Module Info: need xyt.{c,h} gs.{c,h} comm.{c,h} ivec.{c,h} error.{c,h}

347: author:  Henry M. Tufo III
348: e-mail:  hmt@asci.uchicago.edu
349: contact:
350: +--------------------------------+--------------------------------+
351: |MCS Division - Building 221     |Department of Computer Science  |
352: |Argonne National Laboratory     |Ryerson 152                     |
353: |9700 S. Cass Avenue             |The University of Chicago       |
354: |Argonne, IL  60439              |Chicago, IL  60637              |
355: |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
356: +--------------------------------+--------------------------------+

358: Last Modification: 3.20.01
359: **************************************xyt.h***********************************/

361: typedef struct xyt_CDT *xyt_ADT;

363: /*************************************xyt.h************************************
364: Function: XYT_new()

366: Return: ADT ptr or NULL upon failure.
367: Description: This function allocates and returns an xyt handle
368: Usage: xyt_handle = xyt_new();
369: **************************************xyt.h***********************************/
370: PETSC_INTERN xyt_ADT XYT_new(void);

372: /*************************************xyt.h************************************
373: Function: XYT_free()

375: Input : pointer to ADT.
376: Description: This function frees the storage associated with an xyt handle
377: Usage: XYT_free(xyt_handle);
378: **************************************xyt.h***********************************/
379: PETSC_INTERN PetscErrorCode XYT_free(xyt_ADT);

381: /*************************************xyt.h************************************
382: Function: XYT_factor

384: Input : ADT ptr,  and pointer to object
385: Output:
386: Return: 0 on failure, 1 on success
387: Description: This function sets the xyt solver

389: xyt assumptions: given n rows of global coarse matrix (E_loc) where
390:    o global dofs N = sum_p(n), p=0,P-1
391:    (i.e. row dist. with no dof replication)
392:    (5.21.00 will handle dif replication case)
393:    o m is the number of columns in E_loc (m>=n)
394:    o local2global holds global number of column i (i=0,...,m-1)
395:    o local2global holds global number of row    i (i=0,...,n-1)
396:    o mylocmatvec performs E_loc . x_loc where x_loc is an vector of
397:    length m in 1-1 correspondence with local2global
398:    (note that gs package takes care of communication).
399:    (note do not zero out upper m-n entries!)
400:    o mylocmatvec(void *grid_data, double *in, double *out)

402: ML beliefs/usage: move this to ML_XYT_factor routine
403:    o my_ml holds address of ML struct associated w/E_loc, grid_data, grid_tag
404:    o grid_tag, grid_data, my_ml used in
405:       ML_Set_CSolve(my_ml, grid_tag, grid_data, ML_Do_CoarseDirect);
406:    o grid_data used in
407:       A_matvec(grid_data,v,u);

409: Usage:
410: **************************************xyt.h***********************************/
411: PETSC_INTERN PetscErrorCode XYT_factor(xyt_ADT,                                                  /* prev. allocated xyt  handle */
412:                                        PetscInt *,                                               /* global column mapping       */
413:                                        PetscInt,                                                 /* local num rows              */
414:                                        PetscInt,                                                 /* local num cols              */
415:                                        PetscErrorCode (*)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc         */
416:                                        void *);                                                  /* grid data for matvec        */

418: /*************************************xyt.h************************************
419: Function: XYT_solve

421: Input : ADT ptr, b (rhs)
422: Output: x (soln)
423: Return:
424: Description: This function performs x = E^-1.b
425: Usage: XYT_solve(xyt_handle, double *x, double *b)
426: **************************************xyt.h***********************************/
427: PETSC_INTERN PetscErrorCode XYT_solve(xyt_ADT, PetscScalar *, PetscScalar *);

429: /*************************************xyt.h************************************
430: Function: XYT_stats

432: Input : handle
433: **************************************xyt.h***********************************/
434: PETSC_INTERN PetscErrorCode XYT_stats(xyt_ADT);

436: /********************************bit_mask.h************************************

438: Author: Henry M. Tufo III

440: e-mail: hmt@cs.brown.edu

442: snail-mail:
443: Division of Applied Mathematics
444: Brown University
445: Providence, RI 02912

447: Last Modification:
448: 11.21.97
449: *********************************bit_mask.h***********************************/
450: PETSC_INTERN PetscInt       PCTFS_div_ceil(PetscInt, PetscInt);
451: PETSC_INTERN PetscErrorCode PCTFS_set_bit_mask(PetscInt *, PetscInt, PetscInt);
452: PETSC_INTERN PetscInt       PCTFS_len_bit_mask(PetscInt);
453: PETSC_INTERN PetscInt       PCTFS_ct_bits(char *, PetscInt);
454: PETSC_INTERN PetscErrorCode PCTFS_bm_to_proc(char *, PetscInt, PetscInt *);
455: PETSC_INTERN PetscInt       PCTFS_len_buf(PetscInt, PetscInt);