xref: /petsc/src/mat/utils/freespace.c (revision 5bd1e5768a3b141382a3ebb2e780cd2d3b3bfbf0)
1be1d678aSKris Buschelman 
2c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
370f19b1fSKris Buschelman 
4a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceGet(PetscInt n,PetscFreeSpaceList *list)
52e111b49SBarry Smith {
6a1a86e44SBarry Smith   PetscFreeSpaceList a;
7dfbe8321SBarry Smith   PetscErrorCode     ierr;
870f19b1fSKris Buschelman 
970f19b1fSKris Buschelman   PetscFunctionBegin;
1095dccacaSBarry Smith   ierr = PetscNew(&a);CHKERRQ(ierr);
11785e854fSJed Brown   ierr = PetscMalloc1(n,&(a->array_head));CHKERRQ(ierr);
128865f1eaSKarl Rupp 
1370f19b1fSKris Buschelman   a->array            = a->array_head;
142e111b49SBarry Smith   a->local_remaining  = n;
1570f19b1fSKris Buschelman   a->local_used       = 0;
1670f19b1fSKris Buschelman   a->total_array_size = 0;
170298fd71SBarry Smith   a->more_space       = NULL;
1870f19b1fSKris Buschelman 
1970f19b1fSKris Buschelman   if (*list) {
2070f19b1fSKris Buschelman     (*list)->more_space = a;
2170f19b1fSKris Buschelman     a->total_array_size = (*list)->total_array_size;
2270f19b1fSKris Buschelman   }
2370f19b1fSKris Buschelman 
242e111b49SBarry Smith   a->total_array_size += n;
2570f19b1fSKris Buschelman   *list                =  a;
2670f19b1fSKris Buschelman   PetscFunctionReturn(0);
2770f19b1fSKris Buschelman }
2870f19b1fSKris Buschelman 
29a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head,PetscInt *space)
3038baddfdSBarry Smith {
31a1a86e44SBarry Smith   PetscFreeSpaceList a;
32dfbe8321SBarry Smith   PetscErrorCode     ierr;
3370f19b1fSKris Buschelman 
3470f19b1fSKris Buschelman   PetscFunctionBegin;
35c05d87d6SBarry Smith   while ((*head)) {
3670f19b1fSKris Buschelman     a      =  (*head)->more_space;
37580bdb30SBarry Smith     ierr   =  PetscArraycpy(space,(*head)->array_head,(*head)->local_used);CHKERRQ(ierr);
3870f19b1fSKris Buschelman     space += (*head)->local_used;
3970f19b1fSKris Buschelman     ierr   =  PetscFree((*head)->array_head);CHKERRQ(ierr);
4070f19b1fSKris Buschelman     ierr   =  PetscFree(*head);CHKERRQ(ierr);
4170f19b1fSKris Buschelman     *head  =  a;
4270f19b1fSKris Buschelman   }
4370f19b1fSKris Buschelman   PetscFunctionReturn(0);
4470f19b1fSKris Buschelman }
457a48dd6fSHong Zhang 
4630cb48eeSHong Zhang /*
47783ef271SHong Zhang   PetscFreeSpaceContiguous_LU -
48783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array
49783ef271SHong Zhang   that enables an efficient matrix triangular solve.
5030cb48eeSHong Zhang 
5130cb48eeSHong Zhang    Input Parameters:
52783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ILU or LU factorization
53*5bd1e576SStefano Zampini .  space - an allocated array with length nnz of factored matrix.
5430cb48eeSHong Zhang .  n - order of the matrix
552ce24eb6SHong Zhang .  bi - row pointer of factored matrix L with length n+1.
56*5bd1e576SStefano 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.
5730cb48eeSHong Zhang 
5830cb48eeSHong Zhang    Output Parameter:
59*5bd1e576SStefano Zampini .  space - column indices are copied into this array with contiguous layout of L and U
60783ef271SHong Zhang 
612ce24eb6SHong Zhang    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed data structure of L and U
6230cb48eeSHong Zhang */
63783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
6412b5cbf3SHong Zhang {
6512b5cbf3SHong Zhang   PetscFreeSpaceList a;
6612b5cbf3SHong Zhang   PetscErrorCode     ierr;
67f268cba8SShri Abhyankar   PetscInt           row,nnz,*bj,*array,total,bi_temp;
68f268cba8SShri Abhyankar   PetscInt           nnzL,nnzU;
69f268cba8SShri Abhyankar 
70f268cba8SShri Abhyankar   PetscFunctionBegin;
71f268cba8SShri Abhyankar   bi_temp = bi[n];
72f268cba8SShri Abhyankar   row     = 0;
73f268cba8SShri Abhyankar   total   = 0;
74f268cba8SShri Abhyankar   nnzL    = bdiag[0];
75c722e25dSBarry Smith   while ((*head)) {
76f268cba8SShri Abhyankar     total += (*head)->local_used;
77f268cba8SShri Abhyankar     array  = (*head)->array_head;
78f268cba8SShri Abhyankar 
790e97c5f4SHong Zhang     while (row < n) {
800e97c5f4SHong Zhang       if (bi[row+1] > total) break;
81f268cba8SShri Abhyankar       /* copy array entries into bj for this row */
82f268cba8SShri Abhyankar       nnz = bi[row+1] - bi[row];
83f268cba8SShri Abhyankar       /* set bi[row] for new datastruct */
84f268cba8SShri Abhyankar       if (row == 0) {
85f268cba8SShri Abhyankar         bi[row] = 0;
86f268cba8SShri Abhyankar       } else {
87f268cba8SShri Abhyankar         bi[row] = bi[row-1] + nnzL; /* nnzL of previous row */
88f268cba8SShri Abhyankar       }
89f268cba8SShri Abhyankar 
90f268cba8SShri Abhyankar       /* L part */
91f268cba8SShri Abhyankar       nnzL = bdiag[row];
92f268cba8SShri Abhyankar       bj   = space+bi[row];
93580bdb30SBarry Smith       ierr = PetscArraycpy(bj,array,nnzL);CHKERRQ(ierr);
94f268cba8SShri Abhyankar 
95f268cba8SShri Abhyankar       /* diagonal entry */
96f268cba8SShri Abhyankar       bdiag[row]        = bi_temp - 1;
97f268cba8SShri Abhyankar       space[bdiag[row]] = row;
98f268cba8SShri Abhyankar 
99f268cba8SShri Abhyankar       /* U part */
100f268cba8SShri Abhyankar       nnzU    = nnz - nnzL;
101f268cba8SShri Abhyankar       bi_temp = bi_temp - nnzU;
102f268cba8SShri Abhyankar       nnzU--;       /* exclude diagonal */
103f268cba8SShri Abhyankar       bj     = space + bi_temp;
104580bdb30SBarry Smith       ierr   = PetscArraycpy(bj,array+nnzL+1,nnzU);CHKERRQ(ierr);
105f268cba8SShri Abhyankar       array += nnz;
106f268cba8SShri Abhyankar       row++;
107f268cba8SShri Abhyankar     }
108f268cba8SShri Abhyankar 
109f268cba8SShri Abhyankar     a     = (*head)->more_space;
110f268cba8SShri Abhyankar     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
111f268cba8SShri Abhyankar     ierr  = PetscFree(*head);CHKERRQ(ierr);
112f268cba8SShri Abhyankar     *head = a;
113f268cba8SShri Abhyankar   }
11443c97eaeSBarry Smith   if (n) {
115f268cba8SShri Abhyankar     bi[n]    = bi[n-1] + nnzL;
116f268cba8SShri Abhyankar     bdiag[n] = bdiag[n-1]-1;
11743c97eaeSBarry Smith   }
118f268cba8SShri Abhyankar   PetscFunctionReturn(0);
119f268cba8SShri Abhyankar }
120f268cba8SShri Abhyankar 
121783ef271SHong Zhang /*
122783ef271SHong Zhang   PetscFreeSpaceContiguous_Cholesky -
123783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array
124783ef271SHong Zhang   that enables an efficient matrix triangular solve.
125783ef271SHong Zhang 
126783ef271SHong Zhang    Input Parameters:
127783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization
128*5bd1e576SStefano Zampini .  space - an allocated array with length nnz of factored matrix.
129783ef271SHong Zhang .  n - order of the matrix
130783ef271SHong Zhang .  ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix.
131*5bd1e576SStefano Zampini -  udiag - array of length n.
132783ef271SHong Zhang 
133783ef271SHong Zhang    Output Parameter:
134*5bd1e576SStefano Zampini +  space - column indices are copied into this array with contiguous layout of U, with diagonal located as the last entry in each row
135783ef271SHong Zhang -  udiag - indices of diagonal entries
136783ef271SHong Zhang 
137783ef271SHong Zhang    See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description.
138783ef271SHong Zhang */
139783ef271SHong Zhang 
140783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *ui,PetscInt *udiag)
141783ef271SHong Zhang {
142783ef271SHong Zhang   PetscFreeSpaceList a;
143783ef271SHong Zhang   PetscErrorCode     ierr;
144783ef271SHong Zhang   PetscInt           row,nnz,*uj,*array,total;
145783ef271SHong Zhang 
146783ef271SHong Zhang   PetscFunctionBegin;
147783ef271SHong Zhang   row   = 0;
148783ef271SHong Zhang   total = 0;
1490e97c5f4SHong Zhang   while (*head) {
150783ef271SHong Zhang     total += (*head)->local_used;
151783ef271SHong Zhang     array  = (*head)->array_head;
152783ef271SHong Zhang 
1530e97c5f4SHong Zhang     while (row < n) {
1540e97c5f4SHong Zhang       if (ui[row+1] > total) break;
155783ef271SHong Zhang       udiag[row] = ui[row+1] - 1;     /* points to the last entry of U(row,:) */
156783ef271SHong Zhang       nnz        = ui[row+1] - ui[row] - 1; /* exclude diagonal */
157783ef271SHong Zhang       uj         = space + ui[row];
158580bdb30SBarry Smith       ierr       = PetscArraycpy(uj,array+1,nnz);CHKERRQ(ierr);
159783ef271SHong Zhang       uj[nnz]    = array[0]; /* diagonal */
160783ef271SHong Zhang       array     += nnz + 1;
161783ef271SHong Zhang       row++;
162783ef271SHong Zhang     }
163783ef271SHong Zhang 
164783ef271SHong Zhang     a     = (*head)->more_space;
165783ef271SHong Zhang     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
166783ef271SHong Zhang     ierr  = PetscFree(*head);CHKERRQ(ierr);
167783ef271SHong Zhang     *head = a;
168783ef271SHong Zhang   }
169783ef271SHong Zhang   PetscFunctionReturn(0);
170783ef271SHong Zhang }
171783ef271SHong Zhang 
172a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
1737a48dd6fSHong Zhang {
174a1a86e44SBarry Smith   PetscFreeSpaceList a;
1757a48dd6fSHong Zhang   PetscErrorCode     ierr;
1767a48dd6fSHong Zhang 
1777a48dd6fSHong Zhang   PetscFunctionBegin;
178c722e25dSBarry Smith   while ((head)) {
1797a48dd6fSHong Zhang     a    = (head)->more_space;
1807a48dd6fSHong Zhang     ierr = PetscFree((head)->array_head);CHKERRQ(ierr);
1817a48dd6fSHong Zhang     ierr = PetscFree(head);CHKERRQ(ierr);
1827a48dd6fSHong Zhang     head = a;
1837a48dd6fSHong Zhang   }
1847a48dd6fSHong Zhang   PetscFunctionReturn(0);
1857a48dd6fSHong Zhang }
186