xref: /petsc/src/mat/utils/freespace.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1be1d678aSKris Buschelman 
2c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
370f19b1fSKris Buschelman 
4d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceGet(PetscInt n, PetscFreeSpaceList *list)
5d71ae5a4SJacob Faibussowitsch {
6a1a86e44SBarry Smith   PetscFreeSpaceList a;
770f19b1fSKris Buschelman 
870f19b1fSKris Buschelman   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscNew(&a));
109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &(a->array_head)));
118865f1eaSKarl Rupp 
1270f19b1fSKris Buschelman   a->array            = a->array_head;
132e111b49SBarry Smith   a->local_remaining  = n;
1470f19b1fSKris Buschelman   a->local_used       = 0;
1570f19b1fSKris Buschelman   a->total_array_size = 0;
160298fd71SBarry Smith   a->more_space       = NULL;
1770f19b1fSKris Buschelman 
1870f19b1fSKris Buschelman   if (*list) {
1970f19b1fSKris Buschelman     (*list)->more_space = a;
2070f19b1fSKris Buschelman     a->total_array_size = (*list)->total_array_size;
2170f19b1fSKris Buschelman   }
2270f19b1fSKris Buschelman 
232e111b49SBarry Smith   a->total_array_size += n;
2470f19b1fSKris Buschelman   *list = a;
25*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2670f19b1fSKris Buschelman }
2770f19b1fSKris Buschelman 
28d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head, PetscInt *space)
29d71ae5a4SJacob Faibussowitsch {
30a1a86e44SBarry Smith   PetscFreeSpaceList a;
3170f19b1fSKris Buschelman 
3270f19b1fSKris Buschelman   PetscFunctionBegin;
33c05d87d6SBarry Smith   while ((*head)) {
3470f19b1fSKris Buschelman     a = (*head)->more_space;
359566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(space, (*head)->array_head, (*head)->local_used));
3670f19b1fSKris Buschelman     space += (*head)->local_used;
379566063dSJacob Faibussowitsch     PetscCall(PetscFree((*head)->array_head));
389566063dSJacob Faibussowitsch     PetscCall(PetscFree(*head));
3970f19b1fSKris Buschelman     *head = a;
4070f19b1fSKris Buschelman   }
41*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4270f19b1fSKris Buschelman }
437a48dd6fSHong Zhang 
4430cb48eeSHong Zhang /*
45783ef271SHong Zhang   PetscFreeSpaceContiguous_LU -
46783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array
47783ef271SHong Zhang   that enables an efficient matrix triangular solve.
4830cb48eeSHong Zhang 
4930cb48eeSHong Zhang    Input Parameters:
50783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ILU or LU factorization
515bd1e576SStefano Zampini .  space - an allocated array with length nnz of factored matrix.
5230cb48eeSHong Zhang .  n - order of the matrix
532ce24eb6SHong Zhang .  bi - row pointer of factored matrix L with length n+1.
545bd1e576SStefano 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.
5530cb48eeSHong Zhang 
5630cb48eeSHong Zhang    Output Parameter:
575bd1e576SStefano Zampini .  space - column indices are copied into this array with contiguous layout of L and U
58783ef271SHong Zhang 
592ce24eb6SHong Zhang    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed data structure of L and U
6030cb48eeSHong Zhang */
61d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head, PetscInt *space, PetscInt n, PetscInt *bi, PetscInt *bdiag)
62d71ae5a4SJacob Faibussowitsch {
6312b5cbf3SHong Zhang   PetscFreeSpaceList a;
64f268cba8SShri Abhyankar   PetscInt           row, nnz, *bj, *array, total, bi_temp;
65f268cba8SShri Abhyankar   PetscInt           nnzL, nnzU;
66f268cba8SShri Abhyankar 
67f268cba8SShri Abhyankar   PetscFunctionBegin;
68f268cba8SShri Abhyankar   bi_temp = bi[n];
69f268cba8SShri Abhyankar   row     = 0;
70f268cba8SShri Abhyankar   total   = 0;
71f268cba8SShri Abhyankar   nnzL    = bdiag[0];
72c722e25dSBarry Smith   while ((*head)) {
73f268cba8SShri Abhyankar     total += (*head)->local_used;
74f268cba8SShri Abhyankar     array = (*head)->array_head;
75f268cba8SShri Abhyankar 
760e97c5f4SHong Zhang     while (row < n) {
770e97c5f4SHong Zhang       if (bi[row + 1] > total) break;
78f268cba8SShri Abhyankar       /* copy array entries into bj for this row */
79f268cba8SShri Abhyankar       nnz = bi[row + 1] - bi[row];
80f268cba8SShri Abhyankar       /* set bi[row] for new datastruct */
81f268cba8SShri Abhyankar       if (row == 0) {
82f268cba8SShri Abhyankar         bi[row] = 0;
83f268cba8SShri Abhyankar       } else {
84f268cba8SShri Abhyankar         bi[row] = bi[row - 1] + nnzL; /* nnzL of previous row */
85f268cba8SShri Abhyankar       }
86f268cba8SShri Abhyankar 
87f268cba8SShri Abhyankar       /* L part */
88f268cba8SShri Abhyankar       nnzL = bdiag[row];
89f268cba8SShri Abhyankar       bj   = space + bi[row];
909566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(bj, array, nnzL));
91f268cba8SShri Abhyankar 
92f268cba8SShri Abhyankar       /* diagonal entry */
93f268cba8SShri Abhyankar       bdiag[row]        = bi_temp - 1;
94f268cba8SShri Abhyankar       space[bdiag[row]] = row;
95f268cba8SShri Abhyankar 
96f268cba8SShri Abhyankar       /* U part */
97f268cba8SShri Abhyankar       nnzU    = nnz - nnzL;
98f268cba8SShri Abhyankar       bi_temp = bi_temp - nnzU;
99f268cba8SShri Abhyankar       nnzU--; /* exclude diagonal */
100f268cba8SShri Abhyankar       bj = space + bi_temp;
1019566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(bj, array + nnzL + 1, nnzU));
102f268cba8SShri Abhyankar       array += nnz;
103f268cba8SShri Abhyankar       row++;
104f268cba8SShri Abhyankar     }
105f268cba8SShri Abhyankar 
106f268cba8SShri Abhyankar     a = (*head)->more_space;
1079566063dSJacob Faibussowitsch     PetscCall(PetscFree((*head)->array_head));
1089566063dSJacob Faibussowitsch     PetscCall(PetscFree(*head));
109f268cba8SShri Abhyankar     *head = a;
110f268cba8SShri Abhyankar   }
11143c97eaeSBarry Smith   if (n) {
112f268cba8SShri Abhyankar     bi[n]    = bi[n - 1] + nnzL;
113f268cba8SShri Abhyankar     bdiag[n] = bdiag[n - 1] - 1;
11443c97eaeSBarry Smith   }
115*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116f268cba8SShri Abhyankar }
117f268cba8SShri Abhyankar 
118783ef271SHong Zhang /*
119783ef271SHong Zhang   PetscFreeSpaceContiguous_Cholesky -
120783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array
121783ef271SHong Zhang   that enables an efficient matrix triangular solve.
122783ef271SHong Zhang 
123783ef271SHong Zhang    Input Parameters:
124783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization
1255bd1e576SStefano Zampini .  space - an allocated array with length nnz of factored matrix.
126783ef271SHong Zhang .  n - order of the matrix
127783ef271SHong Zhang .  ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix.
1285bd1e576SStefano Zampini -  udiag - array of length n.
129783ef271SHong Zhang 
130783ef271SHong Zhang    Output Parameter:
1315bd1e576SStefano Zampini +  space - column indices are copied into this array with contiguous layout of U, with diagonal located as the last entry in each row
132783ef271SHong Zhang -  udiag - indices of diagonal entries
133783ef271SHong Zhang 
134783ef271SHong Zhang    See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description.
135783ef271SHong Zhang */
136783ef271SHong Zhang 
137d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head, PetscInt *space, PetscInt n, PetscInt *ui, PetscInt *udiag)
138d71ae5a4SJacob Faibussowitsch {
139783ef271SHong Zhang   PetscFreeSpaceList a;
140783ef271SHong Zhang   PetscInt           row, nnz, *uj, *array, total;
141783ef271SHong Zhang 
142783ef271SHong Zhang   PetscFunctionBegin;
143783ef271SHong Zhang   row   = 0;
144783ef271SHong Zhang   total = 0;
1450e97c5f4SHong Zhang   while (*head) {
146783ef271SHong Zhang     total += (*head)->local_used;
147783ef271SHong Zhang     array = (*head)->array_head;
148783ef271SHong Zhang 
1490e97c5f4SHong Zhang     while (row < n) {
1500e97c5f4SHong Zhang       if (ui[row + 1] > total) break;
151783ef271SHong Zhang       udiag[row] = ui[row + 1] - 1;           /* points to the last entry of U(row,:) */
152783ef271SHong Zhang       nnz        = ui[row + 1] - ui[row] - 1; /* exclude diagonal */
153783ef271SHong Zhang       uj         = space + ui[row];
1549566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(uj, array + 1, nnz));
155783ef271SHong Zhang       uj[nnz] = array[0]; /* diagonal */
156783ef271SHong Zhang       array += nnz + 1;
157783ef271SHong Zhang       row++;
158783ef271SHong Zhang     }
159783ef271SHong Zhang 
160783ef271SHong Zhang     a = (*head)->more_space;
1619566063dSJacob Faibussowitsch     PetscCall(PetscFree((*head)->array_head));
1629566063dSJacob Faibussowitsch     PetscCall(PetscFree(*head));
163783ef271SHong Zhang     *head = a;
164783ef271SHong Zhang   }
165*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166783ef271SHong Zhang }
167783ef271SHong Zhang 
168d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
169d71ae5a4SJacob Faibussowitsch {
170a1a86e44SBarry Smith   PetscFreeSpaceList a;
1717a48dd6fSHong Zhang 
1727a48dd6fSHong Zhang   PetscFunctionBegin;
173c722e25dSBarry Smith   while ((head)) {
1747a48dd6fSHong Zhang     a = (head)->more_space;
1759566063dSJacob Faibussowitsch     PetscCall(PetscFree((head)->array_head));
1769566063dSJacob Faibussowitsch     PetscCall(PetscFree(head));
1777a48dd6fSHong Zhang     head = a;
1787a48dd6fSHong Zhang   }
179*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1807a48dd6fSHong Zhang }
181