Actual source code: freespace.c
1: #define PETSCMAT_DLL
3: #include ../src/mat/utils/freespace.h
7: PetscErrorCode PetscFreeSpaceGet(PetscInt n,PetscFreeSpaceList *list)
8: {
9: PetscFreeSpaceList a;
10: PetscErrorCode ierr;
13: PetscMalloc(sizeof(struct _Space),&a);
14: PetscMalloc(n*sizeof(PetscInt),&(a->array_head));
15: a->array = a->array_head;
16: a->local_remaining = n;
17: a->local_used = 0;
18: a->total_array_size = 0;
19: a->more_space = NULL;
21: if (*list) {
22: (*list)->more_space = a;
23: a->total_array_size = (*list)->total_array_size;
24: }
26: a->total_array_size += n;
27: *list = a;
28: return(0);
29: }
33: PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head,PetscInt *space)
34: {
35: PetscFreeSpaceList a;
36: PetscErrorCode ierr;
39: while ((*head)) {
40: a = (*head)->more_space;
41: PetscMemcpy(space,(*head)->array_head,((*head)->local_used)*sizeof(PetscInt));
42: space += (*head)->local_used;
43: PetscFree((*head)->array_head);
44: PetscFree(*head);
45: *head = a;
46: }
47: return(0);
48: }
50: /*
51: PetscFreeSpaceContiguous_LU -
52: Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array
53: that enables an efficient matrix triangular solve.
55: Input Parameters:
56: + head - linked list of column indices obtained from matrix symbolic ILU or LU factorization
57: . space - an allocated int array with length nnz of factored matrix.
58: . n - order of the matrix
59: . bi - row pointer of factored matrix L with length n+1.
60: - bdiag - int array of length n+1. bdiag[i] points to diagonal of U(i,:), and bdiag[n] points to entry of U(n-1,0)-1.
62: Output Parameter:
63: . space - column indices are copied into this int array with contiguous layout of L and U
65: See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed data structure of L and U
66: */
69: PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
70: {
71: PetscFreeSpaceList a;
72: PetscErrorCode ierr;
73: PetscInt row,nnz,*bj,*array,total,bi_temp;
74: PetscInt nnzL,nnzU;
77: bi_temp = bi[n];
78: row = 0;
79: total = 0;
80: nnzL = bdiag[0];
81: while ((*head)!=NULL) {
82: total += (*head)->local_used;
83: array = (*head)->array_head;
84:
85: while (row < n) {
86: if (bi[row+1] > total) break;
87: /* copy array entries into bj for this row */
88: nnz = bi[row+1] - bi[row];
89: /* set bi[row] for new datastruct */
90: if (row == 0 ){
91: bi[row] = 0;
92: } else {
93: bi[row] = bi[row-1] + nnzL; /* nnzL of previous row */
94: }
96: /* L part */
97: nnzL = bdiag[row];
98: bj = space+bi[row];
99: PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));
100:
101: /* diagonal entry */
102: bdiag[row] = bi_temp - 1;
103: space[bdiag[row]] = row;
105: /* U part */
106: nnzU = nnz - nnzL;
107: bi_temp = bi_temp - nnzU;
108: nnzU --; /* exclude diagonal */
109: bj = space + bi_temp;
110: PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));
111: array += nnz;
112: row++;
113: }
115: a = (*head)->more_space;
116: PetscFree((*head)->array_head);
117: PetscFree(*head);
118: *head = a;
119: }
120: bi[n] = bi[n-1] + nnzL;
121: bdiag[n] = bdiag[n-1]-1;
122: return(0);
123: }
125: /*
126: PetscFreeSpaceContiguous_Cholesky -
127: Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array
128: that enables an efficient matrix triangular solve.
130: Input Parameters:
131: + head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization
132: . space - an allocated int array with length nnz of factored matrix.
133: . n - order of the matrix
134: . ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix.
135: - udiag - int array of length n.
137: Output Parameter:
138: + space - column indices are copied into this int array with contiguous layout of U, with diagonal located as the last entry in each row
139: - udiag - indices of diagonal entries
141: See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description.
142: */
146: PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *ui,PetscInt *udiag)
147: {
148: PetscFreeSpaceList a;
149: PetscErrorCode ierr;
150: PetscInt row,nnz,*uj,*array,total;
153: row = 0;
154: total = 0;
155: while (*head) {
156: total += (*head)->local_used;
157: array = (*head)->array_head;
158:
159: while (row < n){
160: if (ui[row+1] > total) break;
161: udiag[row] = ui[row+1] - 1; /* points to the last entry of U(row,:) */
162: nnz = ui[row+1] - ui[row] - 1; /* exclude diagonal */
163: uj = space + ui[row];
164: PetscMemcpy(uj,array+1,nnz*sizeof(PetscInt));
165: uj[nnz] = array[0]; /* diagonal */
166: array += nnz + 1;
167: row++;
168: }
170: a = (*head)->more_space;
171: PetscFree((*head)->array_head);
172: PetscFree(*head);
173: *head = a;
174: }
175: return(0);
176: }
180: PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
181: {
182: PetscFreeSpaceList a;
183: PetscErrorCode ierr;
186: while ((head)!=NULL) {
187: a = (head)->more_space;
188: PetscFree((head)->array_head);
189: PetscFree(head);
190: head = a;
191: }
192: return(0);
193: }