xref: /petsc/include/petscsys.h (revision f68b968ce39302dfa79eb1a6cfa1998ce074e829)
1 /*
2     Provides access to system related and general utility routines.
3 */
4 #if !defined(__PETSCSYS_H)
5 #define __PETSCSYS_H
6 #include "petsc.h"
7 PETSC_EXTERN_CXX_BEGIN
8 
9 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetArchType(char[],size_t);
10 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetHostName(char[],size_t);
11 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetUserName(char[],size_t);
12 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetProgramName(char[],size_t);
13 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSetProgramName(const char[]);
14 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetDate(char[],size_t);
15 
16 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSortInt(PetscInt,PetscInt[]);
17 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
18 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
19 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
20 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
21 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSortReal(PetscInt,PetscReal[]);
22 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
23 
24 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscSetDisplay(void);
25 EXTERN PetscErrorCode PETSC_DLLEXPORT  PetscGetDisplay(char[],size_t);
26 
27 extern PetscCookie PETSC_DLLEXPORT PETSC_RANDOM_COOKIE;
28 
29 typedef enum {RANDOM_DEFAULT,RANDOM_DEFAULT_REAL,RANDOM_DEFAULT_IMAGINARY} PetscRandomType;
30 
31 /*S
32      PetscRandom - Abstract PETSc object that manages generating random numbers
33 
34    Level: intermediate
35 
36   Concepts: random numbers
37 
38 .seealso:  PetscRandomCreate(), PetscRandomGetValue()
39 S*/
40 typedef struct _p_PetscRandom*   PetscRandom;
41 
42 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomCreate(MPI_Comm,PetscRandomType,PetscRandom*);
43 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomGetValue(PetscRandom,PetscScalar*);
44 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
45 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
46 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomSetSeed(PetscRandom,unsigned long);
47 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomGetSeed(PetscRandom,unsigned long *);
48 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomSeed(PetscRandom);
49 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscRandomDestroy(PetscRandom);
50 
51 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetFullPath(const char[],char[],size_t);
52 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetRelativePath(const char[],char[],size_t);
53 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetWorkingDirectory(char[],size_t);
54 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetRealPath(char[],char[]);
55 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetHomeDirectory(char[],size_t);
56 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscTestFile(const char[],char,PetscTruth*);
57 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscTestDirectory(const char[],char,PetscTruth*);
58 
59 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscBinaryRead(int,void*,PetscInt,PetscDataType);
60 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSynchronizedBinaryRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
61 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSynchronizedBinaryWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscTruth);
62 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscTruth);
63 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscBinaryOpen(const char[],PetscFileMode,int *);
64 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscBinaryClose(int);
65 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSharedTmp(MPI_Comm,PetscTruth *);
66 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSharedWorkingDirectory(MPI_Comm,PetscTruth *);
67 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGetTmp(MPI_Comm,char *,size_t);
68 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscFileRetrieve(MPI_Comm,const char *,char *,size_t,PetscTruth*);
69 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscLs(MPI_Comm,const char[],char*,size_t,PetscTruth*);
70 
71 /*
72    In binary files variables are stored using the following lengths,
73   regardless of how they are stored in memory on any one particular
74   machine. Use these rather then sizeof() in computing sizes for
75   PetscBinarySeek().
76 */
77 #define PETSC_BINARY_INT_SIZE    (32/8)
78 #define PETSC_BINARY_FLOAT_SIZE  (32/8)
79 #define PETSC_BINARY_CHAR_SIZE    (8/8)
80 #define PETSC_BINARY_SHORT_SIZE  (16/8)
81 #define PETSC_BINARY_DOUBLE_SIZE (64/8)
82 #define PETSC_BINARY_SCALAR_SIZE sizeof(PetscScalar)
83 
84 /*E
85   PetscBinarySeekType - argument to PetscBinarySeek()
86 
87   Level: advanced
88 
89 .seealso: PetscBinarySeek(), PetscSynchronizedBinarySeek()
90 E*/
91 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
92 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
93 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSynchronizedBinarySeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
94 
95 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSetDebugger(const char[],PetscTruth);
96 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSetDefaultDebugger(void);
97 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSetDebuggerFromString(char*);
98 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscAttachDebugger(void);
99 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscStopForDebugger(void);
100 
101 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGatherNumberOfMessages(MPI_Comm,PetscMPIInt*,PetscMPIInt*,PetscMPIInt*);
102 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,PetscMPIInt*,PetscMPIInt**,PetscMPIInt**);
103 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,PetscMPIInt*,PetscMPIInt*,PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
104 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,PetscMPIInt*,PetscMPIInt*,PetscInt***,MPI_Request**);
105 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,PetscMPIInt*,PetscMPIInt*,PetscScalar***,MPI_Request**);
106 
107 EXTERN PetscErrorCode PETSC_DLLEXPORT PetscSSEIsEnabled(MPI_Comm,PetscTruth *,PetscTruth *);
108 
109 /* Parallel communication routines */
110 /*E
111   InsertMode - Whether entries are inserted or added into vectors or matrices
112 
113   Level: beginner
114 
115 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
116           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
117           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
118 E*/
119 typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES} InsertMode;
120 
121 /*MC
122     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value
123 
124     Level: beginner
125 
126 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
127           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES, INSERT_VALUES,
128           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
129 
130 M*/
131 
132 /*MC
133     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
134                 value into that location
135 
136     Level: beginner
137 
138 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
139           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES, INSERT_VALUES,
140           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
141 
142 M*/
143 
144 /*MC
145     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location
146 
147     Level: beginner
148 
149 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES
150 
151 M*/
152 
153 /*E
154   ScatterMode - Determines the direction of a scatter
155 
156   Level: beginner
157 
158 .seealso: VecScatter, VecScatterBegin(), VecScatterEnd()
159 E*/
160 typedef enum {SCATTER_FORWARD=0, SCATTER_REVERSE=1, SCATTER_FORWARD_LOCAL=2, SCATTER_REVERSE_LOCAL=3, SCATTER_LOCAL=2} ScatterMode;
161 
162 /*MC
163     SCATTER_FORWARD - Scatters the values as dictated by the VecScatterCreate() call
164 
165     Level: beginner
166 
167 .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_REVERSE, SCATTER_FORWARD_LOCAL,
168           SCATTER_REVERSE_LOCAL
169 
170 M*/
171 
172 /*MC
173     SCATTER_REVERSE - Moves the values in the opposite direction then the directions indicated in
174          in the VecScatterCreate()
175 
176     Level: beginner
177 
178 .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_FORWARD, SCATTER_FORWARD_LOCAL,
179           SCATTER_REVERSE_LOCAL
180 
181 M*/
182 
183 /*MC
184     SCATTER_FORWARD_LOCAL - Scatters the values as dictated by the VecScatterCreate() call except NO parallel communication
185        is done. Any variables that have be moved between processes are ignored
186 
187     Level: developer
188 
189 .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_REVERSE, SCATTER_FORWARD,
190           SCATTER_REVERSE_LOCAL
191 
192 M*/
193 
194 /*MC
195     SCATTER_REVERSE_LOCAL - Moves the values in the opposite direction then the directions indicated in
196          in the VecScatterCreate()  except NO parallel communication
197        is done. Any variables that have be moved between processes are ignored
198 
199     Level: developer
200 
201 .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_FORWARD, SCATTER_FORWARD_LOCAL,
202           SCATTER_REVERSE
203 
204 M*/
205 
206 /*
207   Create and initialize a linked list
208   Input Parameters:
209     idx_start - starting index of the list
210     lnk_max   - max value of lnk indicating the end of the list
211     nlnk      - max length of the list
212   Output Parameters:
213     lnk       - list initialized
214     bt        - PetscBT (bitarray) with all bits set to false
215 */
216 #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
217   (PetscMalloc(nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,bt) || PetscBTMemzero(nlnk,bt) || (lnk[idx_start] = lnk_max,0))
218 
219 /*
220   Add a index set into a sorted linked list
221   Input Parameters:
222     nidx      - number of input indices
223     indices   - interger array
224     idx_start - starting index of the list
225     lnk       - linked list(an integer array) that is created
226     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
227   output Parameters:
228     nlnk      - number of newly added indices
229     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
230     bt        - updated PetscBT (bitarray)
231 */
232 #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
233 {\
234   PetscInt _k,_entry,_location,_lnkdata;\
235   nlnk     = 0;\
236   _lnkdata = idx_start;\
237   for (_k=0; _k<nidx; _k++){\
238     _entry = indices[_k];\
239     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
240       /* search for insertion location */\
241       /* start from the beginning if _entry < previous _entry */\
242       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
243       do {\
244         _location = _lnkdata;\
245         _lnkdata  = lnk[_location];\
246       } while (_entry > _lnkdata);\
247       /* insertion location is found, add entry into lnk */\
248       lnk[_location] = _entry;\
249       lnk[_entry]    = _lnkdata;\
250       nlnk++;\
251       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
252     }\
253   }\
254 }
255 
256 /*
257   Add a permumted index set into a sorted linked list
258   Input Parameters:
259     nidx      - number of input indices
260     indices   - interger array
261     perm      - permutation of indices
262     idx_start - starting index of the list
263     lnk       - linked list(an integer array) that is created
264     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
265   output Parameters:
266     nlnk      - number of newly added indices
267     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
268     bt        - updated PetscBT (bitarray)
269 */
270 #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
271 {\
272   PetscInt _k,_entry,_location,_lnkdata;\
273   nlnk     = 0;\
274   _lnkdata = idx_start;\
275   for (_k=0; _k<nidx; _k++){\
276     _entry = perm[indices[_k]];\
277     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
278       /* search for insertion location */\
279       /* start from the beginning if _entry < previous _entry */\
280       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
281       do {\
282         _location = _lnkdata;\
283         _lnkdata  = lnk[_location];\
284       } while (_entry > _lnkdata);\
285       /* insertion location is found, add entry into lnk */\
286       lnk[_location] = _entry;\
287       lnk[_entry]    = _lnkdata;\
288       nlnk++;\
289       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
290     }\
291   }\
292 }
293 
294 /*
295   Add a SORTED index set into a sorted linked list
296   Input Parameters:
297     nidx      - number of input indices
298     indices   - sorted interger array
299     idx_start - starting index of the list
300     lnk       - linked list(an integer array) that is created
301     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
302   output Parameters:
303     nlnk      - number of newly added indices
304     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
305     bt        - updated PetscBT (bitarray)
306 */
307 #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
308 {\
309   PetscInt _k,_entry,_location,_lnkdata;\
310   nlnk      = 0;\
311   _lnkdata  = idx_start;\
312   for (_k=0; _k<nidx; _k++){\
313     _entry = indices[_k];\
314     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
315       /* search for insertion location */\
316       do {\
317         _location = _lnkdata;\
318         _lnkdata  = lnk[_location];\
319       } while (_entry > _lnkdata);\
320       /* insertion location is found, add entry into lnk */\
321       lnk[_location] = _entry;\
322       lnk[_entry]    = _lnkdata;\
323       nlnk++;\
324       _lnkdata = _entry; /* next search starts from here */\
325     }\
326   }\
327 }
328 
329 /*
330   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
331   Same as PetscLLAddSorted() with an additional operation:
332        count the number of input indices that are no larger than 'diag'
333   Input Parameters:
334     indices   - sorted interger array
335     idx_start - starting index of the list
336     lnk       - linked list(an integer array) that is created
337     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
338     diag      - index of the active row in LUFactorSymbolic
339     nzbd      - number of input indices with indices <= idx_start
340   output Parameters:
341     nlnk      - number of newly added indices
342     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
343     bt        - updated PetscBT (bitarray)
344     im        - im[idx_start] =  num of entries with indices <= diag
345 */
346 #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
347 {\
348   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
349   nlnk     = 0;\
350   _lnkdata = idx_start;\
351   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
352   for (_k=0; _k<_nidx; _k++){\
353     _entry = indices[_k];\
354     nzbd++;\
355     if ( _entry== diag) im[idx_start] = nzbd;\
356     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
357       /* search for insertion location */\
358       do {\
359         _location = _lnkdata;\
360         _lnkdata  = lnk[_location];\
361       } while (_entry > _lnkdata);\
362       /* insertion location is found, add entry into lnk */\
363       lnk[_location] = _entry;\
364       lnk[_entry]    = _lnkdata;\
365       nlnk++;\
366       _lnkdata = _entry; /* next search starts from here */\
367     }\
368   }\
369 }
370 
371 /*
372   Copy data on the list into an array, then initialize the list
373   Input Parameters:
374     idx_start - starting index of the list
375     lnk_max   - max value of lnk indicating the end of the list
376     nlnk      - number of data on the list to be copied
377     lnk       - linked list
378     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
379   output Parameters:
380     indices   - array that contains the copied data
381     lnk       - linked list that is cleaned and initialize
382     bt        - PetscBT (bitarray) with all bits set to false
383 */
384 #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
385 {\
386   PetscInt _j,_idx=idx_start;\
387   for (_j=0; _j<nlnk; _j++){\
388     _idx = lnk[_idx];\
389     *(indices+_j) = _idx;\
390     ierr = PetscBTClear(bt,_idx);CHKERRQ(ierr);\
391   }\
392   lnk[idx_start] = lnk_max;\
393 }
394 /*
395   Free memories used by the list
396 */
397 #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(bt))
398 
399 /* Routines below are used for incomplete matrix factorization */
400 /*
401   Create and initialize a linked list and its levels
402   Input Parameters:
403     idx_start - starting index of the list
404     lnk_max   - max value of lnk indicating the end of the list
405     nlnk      - max length of the list
406   Output Parameters:
407     lnk       - list initialized
408     lnk_lvl   - array of size nlnk for storing levels of lnk
409     bt        - PetscBT (bitarray) with all bits set to false
410 */
411 #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
412   (PetscMalloc(2*nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,bt) || PetscBTMemzero(nlnk,bt) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))
413 
414 /*
415   Initialize a sorted linked list used for ILU and ICC
416   Input Parameters:
417     nidx      - number of input idx
418     idx       - interger array used for storing column indices
419     idx_start - starting index of the list
420     perm      - indices of an IS
421     lnk       - linked list(an integer array) that is created
422     lnklvl    - levels of lnk
423     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
424   output Parameters:
425     nlnk     - number of newly added idx
426     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
427     lnklvl   - levels of lnk
428     bt       - updated PetscBT (bitarray)
429 */
430 #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
431 {\
432   PetscInt _k,_entry,_location,_lnkdata;\
433   nlnk     = 0;\
434   _lnkdata = idx_start;\
435   for (_k=0; _k<nidx; _k++){\
436     _entry = perm[idx[_k]];\
437     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
438       /* search for insertion location */\
439       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
440       do {\
441         _location = _lnkdata;\
442         _lnkdata  = lnk[_location];\
443       } while (_entry > _lnkdata);\
444       /* insertion location is found, add entry into lnk */\
445       lnk[_location]  = _entry;\
446       lnk[_entry]     = _lnkdata;\
447       lnklvl[_entry] = 0;\
448       nlnk++;\
449       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
450     }\
451   }\
452 }
453 
454 /*
455   Add a SORTED index set into a sorted linked list for ILU
456   Input Parameters:
457     nidx      - number of input indices
458     idx       - sorted interger array used for storing column indices
459     level     - level of fill, e.g., ICC(level)
460     idxlvl    - level of idx
461     idx_start - starting index of the list
462     lnk       - linked list(an integer array) that is created
463     lnklvl    - levels of lnk
464     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
465     prow      - the row number of idx
466   output Parameters:
467     nlnk     - number of newly added idx
468     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
469     lnklvl   - levels of lnk
470     bt       - updated PetscBT (bitarray)
471 
472   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
473         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
474 */
475 #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
476 {\
477   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
478   nlnk     = 0;\
479   _lnkdata = idx_start;\
480   for (_k=0; _k<nidx; _k++){\
481     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
482     if (_incrlev > level) continue;\
483     _entry = idx[_k];\
484     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
485       /* search for insertion location */\
486       do {\
487         _location = _lnkdata;\
488         _lnkdata  = lnk[_location];\
489       } while (_entry > _lnkdata);\
490       /* insertion location is found, add entry into lnk */\
491       lnk[_location]  = _entry;\
492       lnk[_entry]     = _lnkdata;\
493       lnklvl[_entry] = _incrlev;\
494       nlnk++;\
495       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
496     } else { /* existing entry: update lnklvl */\
497       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
498     }\
499   }\
500 }
501 
502 /*
503   Add a index set into a sorted linked list
504   Input Parameters:
505     nidx      - number of input idx
506     idx   - interger array used for storing column indices
507     level     - level of fill, e.g., ICC(level)
508     idxlvl - level of idx
509     idx_start - starting index of the list
510     lnk       - linked list(an integer array) that is created
511     lnklvl   - levels of lnk
512     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
513   output Parameters:
514     nlnk      - number of newly added idx
515     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
516     lnklvl   - levels of lnk
517     bt        - updated PetscBT (bitarray)
518 */
519 #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
520 {\
521   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
522   nlnk     = 0;\
523   _lnkdata = idx_start;\
524   for (_k=0; _k<nidx; _k++){\
525     _incrlev = idxlvl[_k] + 1;\
526     if (_incrlev > level) continue;\
527     _entry = idx[_k];\
528     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
529       /* search for insertion location */\
530       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
531       do {\
532         _location = _lnkdata;\
533         _lnkdata  = lnk[_location];\
534       } while (_entry > _lnkdata);\
535       /* insertion location is found, add entry into lnk */\
536       lnk[_location]  = _entry;\
537       lnk[_entry]     = _lnkdata;\
538       lnklvl[_entry] = _incrlev;\
539       nlnk++;\
540       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
541     } else { /* existing entry: update lnklvl */\
542       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
543     }\
544   }\
545 }
546 
547 /*
548   Add a SORTED index set into a sorted linked list
549   Input Parameters:
550     nidx      - number of input indices
551     idx   - sorted interger array used for storing column indices
552     level     - level of fill, e.g., ICC(level)
553     idxlvl - level of idx
554     idx_start - starting index of the list
555     lnk       - linked list(an integer array) that is created
556     lnklvl    - levels of lnk
557     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
558   output Parameters:
559     nlnk      - number of newly added idx
560     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
561     lnklvl    - levels of lnk
562     bt        - updated PetscBT (bitarray)
563 */
564 #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
565 {\
566   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
567   nlnk = 0;\
568   _lnkdata = idx_start;\
569   for (_k=0; _k<nidx; _k++){\
570     _incrlev = idxlvl[_k] + 1;\
571     if (_incrlev > level) continue;\
572     _entry = idx[_k];\
573     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
574       /* search for insertion location */\
575       do {\
576         _location = _lnkdata;\
577         _lnkdata  = lnk[_location];\
578       } while (_entry > _lnkdata);\
579       /* insertion location is found, add entry into lnk */\
580       lnk[_location] = _entry;\
581       lnk[_entry]    = _lnkdata;\
582       lnklvl[_entry] = _incrlev;\
583       nlnk++;\
584       _lnkdata = _entry; /* next search starts from here */\
585     } else { /* existing entry: update lnklvl */\
586       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
587     }\
588   }\
589 }
590 
591 /*
592   Add a SORTED index set into a sorted linked list for ICC
593   Input Parameters:
594     nidx      - number of input indices
595     idx       - sorted interger array used for storing column indices
596     level     - level of fill, e.g., ICC(level)
597     idxlvl    - level of idx
598     idx_start - starting index of the list
599     lnk       - linked list(an integer array) that is created
600     lnklvl    - levels of lnk
601     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
602     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
603   output Parameters:
604     nlnk   - number of newly added indices
605     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
606     lnklvl - levels of lnk
607     bt     - updated PetscBT (bitarray)
608   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
609         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
610 */
611 #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
612 {\
613   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
614   nlnk = 0;\
615   _lnkdata = idx_start;\
616   for (_k=0; _k<nidx; _k++){\
617     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
618     if (_incrlev > level) continue;\
619     _entry = idx[_k];\
620     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
621       /* search for insertion location */\
622       do {\
623         _location = _lnkdata;\
624         _lnkdata  = lnk[_location];\
625       } while (_entry > _lnkdata);\
626       /* insertion location is found, add entry into lnk */\
627       lnk[_location] = _entry;\
628       lnk[_entry]    = _lnkdata;\
629       lnklvl[_entry] = _incrlev;\
630       nlnk++;\
631       _lnkdata = _entry; /* next search starts from here */\
632     } else { /* existing entry: update lnklvl */\
633       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
634     }\
635   }\
636 }
637 
638 /*
639   Copy data on the list into an array, then initialize the list
640   Input Parameters:
641     idx_start - starting index of the list
642     lnk_max   - max value of lnk indicating the end of the list
643     nlnk      - number of data on the list to be copied
644     lnk       - linked list
645     lnklvl    - level of lnk
646     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
647   output Parameters:
648     indices - array that contains the copied data
649     lnk     - linked list that is cleaned and initialize
650     lnklvl  - level of lnk that is reinitialized
651     bt      - PetscBT (bitarray) with all bits set to false
652 */
653 #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
654 {\
655   PetscInt _j,_idx=idx_start;\
656   for (_j=0; _j<nlnk; _j++){\
657     _idx = lnk[_idx];\
658     *(indices+_j) = _idx;\
659     *(indiceslvl+_j) = lnklvl[_idx];\
660     lnklvl[_idx] = -1;\
661     ierr = PetscBTClear(bt,_idx);CHKERRQ(ierr);\
662   }\
663   lnk[idx_start] = lnk_max;\
664 }
665 /*
666   Free memories used by the list
667 */
668 #define PetscIncompleteLLDestroy(lnk,bt) \
669   (PetscFree(lnk) || PetscBTDestroy(bt) || (0))
670 
671 PETSC_EXTERN_CXX_END
672 #endif /* __PETSCSYS_H */
673