xref: /petsc/src/mat/utils/freespace.c (revision 970231d20df44f79b27787157e39d441e79f434b)
1c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
270f19b1fSKris Buschelman 
PetscFreeSpaceGet(PetscInt n,PetscFreeSpaceList * list)3d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceGet(PetscInt n, PetscFreeSpaceList *list)
4d71ae5a4SJacob Faibussowitsch {
5a1a86e44SBarry Smith   PetscFreeSpaceList a;
670f19b1fSKris Buschelman 
770f19b1fSKris Buschelman   PetscFunctionBegin;
89566063dSJacob Faibussowitsch   PetscCall(PetscNew(&a));
9*f4f49eeaSPierre Jolivet   PetscCall(PetscMalloc1(n, &a->array_head));
108865f1eaSKarl Rupp 
1170f19b1fSKris Buschelman   a->array            = a->array_head;
122e111b49SBarry Smith   a->local_remaining  = n;
1370f19b1fSKris Buschelman   a->local_used       = 0;
1470f19b1fSKris Buschelman   a->total_array_size = 0;
150298fd71SBarry Smith   a->more_space       = NULL;
1670f19b1fSKris Buschelman 
1770f19b1fSKris Buschelman   if (*list) {
1870f19b1fSKris Buschelman     (*list)->more_space = a;
1970f19b1fSKris Buschelman     a->total_array_size = (*list)->total_array_size;
2070f19b1fSKris Buschelman   }
2170f19b1fSKris Buschelman 
222e111b49SBarry Smith   a->total_array_size += n;
2370f19b1fSKris Buschelman   *list = a;
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2570f19b1fSKris Buschelman }
2670f19b1fSKris Buschelman 
PetscFreeSpaceContiguous(PetscFreeSpaceList * head,PetscInt * space)27d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head, PetscInt *space)
28d71ae5a4SJacob Faibussowitsch {
29a1a86e44SBarry Smith   PetscFreeSpaceList a;
3070f19b1fSKris Buschelman 
3170f19b1fSKris Buschelman   PetscFunctionBegin;
32*f4f49eeaSPierre Jolivet   while (*head) {
3370f19b1fSKris Buschelman     a = (*head)->more_space;
349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(space, (*head)->array_head, (*head)->local_used));
358e3a54c0SPierre Jolivet     space = PetscSafePointerPlusOffset(space, (*head)->local_used);
369566063dSJacob Faibussowitsch     PetscCall(PetscFree((*head)->array_head));
379566063dSJacob Faibussowitsch     PetscCall(PetscFree(*head));
3870f19b1fSKris Buschelman     *head = a;
3970f19b1fSKris Buschelman   }
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4170f19b1fSKris Buschelman }
427a48dd6fSHong Zhang 
4330cb48eeSHong Zhang /*
44783ef271SHong Zhang   PetscFreeSpaceContiguous_LU -
45783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array
46783ef271SHong Zhang   that enables an efficient matrix triangular solve.
4730cb48eeSHong Zhang 
4830cb48eeSHong Zhang    Input Parameters:
49783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ILU or LU factorization
505bd1e576SStefano Zampini .  space - an allocated array with length nnz of factored matrix.
5130cb48eeSHong Zhang .  n - order of the matrix
522ce24eb6SHong Zhang .  bi - row pointer of factored matrix L with length n+1.
535bd1e576SStefano Zampini -  bdiag - 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.
5430cb48eeSHong Zhang 
5530cb48eeSHong Zhang    Output Parameter:
565bd1e576SStefano Zampini .  space - column indices are copied into this array with contiguous layout of L and U
57783ef271SHong Zhang 
582ce24eb6SHong Zhang    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed data structure of L and U
5930cb48eeSHong Zhang */
PetscFreeSpaceContiguous_LU(PetscFreeSpaceList * head,PetscInt * space,PetscInt n,PetscInt * bi,PetscInt * bdiag)60d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head, PetscInt *space, PetscInt n, PetscInt *bi, PetscInt *bdiag)
61d71ae5a4SJacob Faibussowitsch {
6212b5cbf3SHong Zhang   PetscFreeSpaceList a;
63f268cba8SShri Abhyankar   PetscInt           row, nnz, *bj, *array, total, bi_temp;
64f268cba8SShri Abhyankar   PetscInt           nnzL, nnzU;
65f268cba8SShri Abhyankar 
66f268cba8SShri Abhyankar   PetscFunctionBegin;
67f268cba8SShri Abhyankar   bi_temp = bi[n];
68f268cba8SShri Abhyankar   row     = 0;
69f268cba8SShri Abhyankar   total   = 0;
70f268cba8SShri Abhyankar   nnzL    = bdiag[0];
71*f4f49eeaSPierre Jolivet   while (*head) {
72f268cba8SShri Abhyankar     total += (*head)->local_used;
73f268cba8SShri Abhyankar     array = (*head)->array_head;
74f268cba8SShri Abhyankar 
750e97c5f4SHong Zhang     while (row < n) {
760e97c5f4SHong Zhang       if (bi[row + 1] > total) break;
77f268cba8SShri Abhyankar       /* copy array entries into bj for this row */
78f268cba8SShri Abhyankar       nnz = bi[row + 1] - bi[row];
79f268cba8SShri Abhyankar       /* set bi[row] for new datastruct */
80f268cba8SShri Abhyankar       if (row == 0) {
81f268cba8SShri Abhyankar         bi[row] = 0;
82f268cba8SShri Abhyankar       } else {
83f268cba8SShri Abhyankar         bi[row] = bi[row - 1] + nnzL; /* nnzL of previous row */
84f268cba8SShri Abhyankar       }
85f268cba8SShri Abhyankar 
86f268cba8SShri Abhyankar       /* L part */
87f268cba8SShri Abhyankar       nnzL = bdiag[row];
88f268cba8SShri Abhyankar       bj   = space + bi[row];
899566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(bj, array, nnzL));
90f268cba8SShri Abhyankar 
91f268cba8SShri Abhyankar       /* diagonal entry */
92f268cba8SShri Abhyankar       bdiag[row]        = bi_temp - 1;
93f268cba8SShri Abhyankar       space[bdiag[row]] = row;
94f268cba8SShri Abhyankar 
95f268cba8SShri Abhyankar       /* U part */
96f268cba8SShri Abhyankar       nnzU    = nnz - nnzL;
97f268cba8SShri Abhyankar       bi_temp = bi_temp - nnzU;
98f268cba8SShri Abhyankar       nnzU--; /* exclude diagonal */
99f268cba8SShri Abhyankar       bj = space + bi_temp;
1009566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(bj, array + nnzL + 1, nnzU));
101f268cba8SShri Abhyankar       array += nnz;
102f268cba8SShri Abhyankar       row++;
103f268cba8SShri Abhyankar     }
104f268cba8SShri Abhyankar 
105f268cba8SShri Abhyankar     a = (*head)->more_space;
1069566063dSJacob Faibussowitsch     PetscCall(PetscFree((*head)->array_head));
1079566063dSJacob Faibussowitsch     PetscCall(PetscFree(*head));
108f268cba8SShri Abhyankar     *head = a;
109f268cba8SShri Abhyankar   }
11043c97eaeSBarry Smith   if (n) {
111f268cba8SShri Abhyankar     bi[n]    = bi[n - 1] + nnzL;
112f268cba8SShri Abhyankar     bdiag[n] = bdiag[n - 1] - 1;
11343c97eaeSBarry Smith   }
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115f268cba8SShri Abhyankar }
116f268cba8SShri Abhyankar 
117783ef271SHong Zhang /*
118783ef271SHong Zhang   PetscFreeSpaceContiguous_Cholesky -
119783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array
120783ef271SHong Zhang   that enables an efficient matrix triangular solve.
121783ef271SHong Zhang 
122783ef271SHong Zhang    Input Parameters:
123783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization
1245bd1e576SStefano Zampini .  space - an allocated array with length nnz of factored matrix.
125783ef271SHong Zhang .  n - order of the matrix
126783ef271SHong Zhang .  ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix.
1275bd1e576SStefano Zampini -  udiag - array of length n.
128783ef271SHong Zhang 
129783ef271SHong Zhang    Output Parameter:
1305bd1e576SStefano Zampini +  space - column indices are copied into this array with contiguous layout of U, with diagonal located as the last entry in each row
131783ef271SHong Zhang -  udiag - indices of diagonal entries
132783ef271SHong Zhang 
133783ef271SHong Zhang    See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description.
134783ef271SHong Zhang */
135783ef271SHong Zhang 
PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList * head,PetscInt * space,PetscInt n,PetscInt * ui,PetscInt * udiag)136d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head, PetscInt *space, PetscInt n, PetscInt *ui, PetscInt *udiag)
137d71ae5a4SJacob Faibussowitsch {
138783ef271SHong Zhang   PetscFreeSpaceList a;
139783ef271SHong Zhang   PetscInt           row, nnz, *uj, *array, total;
140783ef271SHong Zhang 
141783ef271SHong Zhang   PetscFunctionBegin;
142783ef271SHong Zhang   row   = 0;
143783ef271SHong Zhang   total = 0;
1440e97c5f4SHong Zhang   while (*head) {
145783ef271SHong Zhang     total += (*head)->local_used;
146783ef271SHong Zhang     array = (*head)->array_head;
147783ef271SHong Zhang 
1480e97c5f4SHong Zhang     while (row < n) {
1490e97c5f4SHong Zhang       if (ui[row + 1] > total) break;
150783ef271SHong Zhang       udiag[row] = ui[row + 1] - 1;           /* points to the last entry of U(row,:) */
151783ef271SHong Zhang       nnz        = ui[row + 1] - ui[row] - 1; /* exclude diagonal */
152783ef271SHong Zhang       uj         = space + ui[row];
1539566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(uj, array + 1, nnz));
154783ef271SHong Zhang       uj[nnz] = array[0]; /* diagonal */
155783ef271SHong Zhang       array += nnz + 1;
156783ef271SHong Zhang       row++;
157783ef271SHong Zhang     }
158783ef271SHong Zhang 
159783ef271SHong Zhang     a = (*head)->more_space;
1609566063dSJacob Faibussowitsch     PetscCall(PetscFree((*head)->array_head));
1619566063dSJacob Faibussowitsch     PetscCall(PetscFree(*head));
162783ef271SHong Zhang     *head = a;
163783ef271SHong Zhang   }
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165783ef271SHong Zhang }
166783ef271SHong Zhang 
PetscFreeSpaceDestroy(PetscFreeSpaceList head)167d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
168d71ae5a4SJacob Faibussowitsch {
169a1a86e44SBarry Smith   PetscFreeSpaceList a;
1707a48dd6fSHong Zhang 
1717a48dd6fSHong Zhang   PetscFunctionBegin;
172*f4f49eeaSPierre Jolivet   while (head) {
1737a48dd6fSHong Zhang     a = (head)->more_space;
1749566063dSJacob Faibussowitsch     PetscCall(PetscFree((head)->array_head));
1759566063dSJacob Faibussowitsch     PetscCall(PetscFree(head));
1767a48dd6fSHong Zhang     head = a;
1777a48dd6fSHong Zhang   }
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1797a48dd6fSHong Zhang }
180