xref: /petsc/src/mat/utils/freespace.c (revision 2ce24eb6ea71247353a331e1f6316651fa3aa6f1)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
37c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
470f19b1fSKris Buschelman 
570f19b1fSKris Buschelman #undef __FUNCT__
6a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceGet"
7a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceGet(PetscInt n,PetscFreeSpaceList *list)
82e111b49SBarry Smith {
9a1a86e44SBarry Smith   PetscFreeSpaceList a;
10dfbe8321SBarry Smith   PetscErrorCode     ierr;
1170f19b1fSKris Buschelman 
1270f19b1fSKris Buschelman   PetscFunctionBegin;
13a1a86e44SBarry Smith   ierr = PetscMalloc(sizeof(struct _Space),&a);CHKERRQ(ierr);
142e111b49SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&(a->array_head));CHKERRQ(ierr);
1570f19b1fSKris Buschelman   a->array            = a->array_head;
162e111b49SBarry Smith   a->local_remaining  = n;
1770f19b1fSKris Buschelman   a->local_used       = 0;
1870f19b1fSKris Buschelman   a->total_array_size = 0;
1970f19b1fSKris Buschelman   a->more_space       = NULL;
2070f19b1fSKris Buschelman 
2170f19b1fSKris Buschelman   if (*list) {
2270f19b1fSKris Buschelman     (*list)->more_space = a;
2370f19b1fSKris Buschelman     a->total_array_size = (*list)->total_array_size;
2470f19b1fSKris Buschelman   }
2570f19b1fSKris Buschelman 
262e111b49SBarry Smith   a->total_array_size += n;
2770f19b1fSKris Buschelman   *list               =  a;
2870f19b1fSKris Buschelman   PetscFunctionReturn(0);
2970f19b1fSKris Buschelman }
3070f19b1fSKris Buschelman 
3170f19b1fSKris Buschelman #undef __FUNCT__
32a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceContiguous"
33a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head,PetscInt *space)
3438baddfdSBarry Smith {
35a1a86e44SBarry Smith   PetscFreeSpaceList a;
36dfbe8321SBarry Smith   PetscErrorCode     ierr;
3770f19b1fSKris Buschelman 
3870f19b1fSKris Buschelman   PetscFunctionBegin;
39c05d87d6SBarry Smith   while ((*head)) {
4070f19b1fSKris Buschelman     a     =  (*head)->more_space;
412e111b49SBarry Smith     ierr  =  PetscMemcpy(space,(*head)->array_head,((*head)->local_used)*sizeof(PetscInt));CHKERRQ(ierr);
4270f19b1fSKris Buschelman     space += (*head)->local_used;
4370f19b1fSKris Buschelman     ierr  =  PetscFree((*head)->array_head);CHKERRQ(ierr);
4470f19b1fSKris Buschelman     ierr  =  PetscFree(*head);CHKERRQ(ierr);
4570f19b1fSKris Buschelman     *head =  a;
4670f19b1fSKris Buschelman   }
4770f19b1fSKris Buschelman   PetscFunctionReturn(0);
4870f19b1fSKris Buschelman }
497a48dd6fSHong Zhang 
5030cb48eeSHong Zhang /*
51783ef271SHong Zhang   PetscFreeSpaceContiguous_LU -
52783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array
53783ef271SHong Zhang   that enables an efficient matrix triangular solve.
5430cb48eeSHong Zhang 
5530cb48eeSHong Zhang    Input Parameters:
56783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ILU or LU factorization
57783ef271SHong Zhang .  space - an allocated int array with length nnz of factored matrix.
5830cb48eeSHong Zhang .  n - order of the matrix
59*2ce24eb6SHong Zhang .  bi - row pointer of factored matrix L with length n+1.
60*2ce24eb6SHong Zhang -  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.
6130cb48eeSHong Zhang 
6230cb48eeSHong Zhang    Output Parameter:
6330cb48eeSHong Zhang .  space - column indices are copied into this int array with contiguous layout of L and U
64783ef271SHong Zhang 
65*2ce24eb6SHong Zhang    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed data structure of L and U
6630cb48eeSHong Zhang */
677a48dd6fSHong Zhang #undef __FUNCT__
68783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_LU"
69783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
7012b5cbf3SHong Zhang {
7112b5cbf3SHong Zhang   PetscFreeSpaceList a;
7212b5cbf3SHong Zhang   PetscErrorCode     ierr;
73f268cba8SShri Abhyankar   PetscInt           row,nnz,*bj,*array,total,bi_temp;
74f268cba8SShri Abhyankar   PetscInt           nnzL,nnzU;
75f268cba8SShri Abhyankar 
76f268cba8SShri Abhyankar   PetscFunctionBegin;
77f268cba8SShri Abhyankar   bi_temp = bi[n];
78f268cba8SShri Abhyankar   row       = 0;
79f268cba8SShri Abhyankar   total     = 0;
80f268cba8SShri Abhyankar   nnzL  = bdiag[0];
81f268cba8SShri Abhyankar   while ((*head)!=NULL) {
82f268cba8SShri Abhyankar     total += (*head)->local_used;
83f268cba8SShri Abhyankar     array  = (*head)->array_head;
84f268cba8SShri Abhyankar 
8581500f6aSJed Brown     while (row < n && bi[row+1] <= total) {
86f268cba8SShri Abhyankar       /* copy array entries into bj for this row */
87f268cba8SShri Abhyankar       nnz  = bi[row+1] - bi[row];
88f268cba8SShri Abhyankar       /* set bi[row] for new datastruct */
89f268cba8SShri Abhyankar       if (row == 0 ){
90f268cba8SShri Abhyankar         bi[row] = 0;
91f268cba8SShri Abhyankar       } else {
92f268cba8SShri Abhyankar         bi[row] = bi[row-1] + nnzL; /* nnzL of previous row */
93f268cba8SShri Abhyankar       }
94f268cba8SShri Abhyankar 
95f268cba8SShri Abhyankar       /* L part */
96f268cba8SShri Abhyankar       nnzL = bdiag[row];
97f268cba8SShri Abhyankar       bj   = space+bi[row];
98f268cba8SShri Abhyankar       ierr = PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));CHKERRQ(ierr);
99f268cba8SShri Abhyankar 
100f268cba8SShri Abhyankar       /* diagonal entry */
101f268cba8SShri Abhyankar       bdiag[row] = bi_temp - 1;
102f268cba8SShri Abhyankar       space[bdiag[row]] = row;
103f268cba8SShri Abhyankar 
104f268cba8SShri Abhyankar       /* U part */
105f268cba8SShri Abhyankar       nnzU        = nnz - nnzL;
106f268cba8SShri Abhyankar       bi_temp = bi_temp - nnzU;
107f268cba8SShri Abhyankar       nnzU --;      /* exclude diagonal */
108f268cba8SShri Abhyankar       bj = space + bi_temp;
109f268cba8SShri Abhyankar       ierr = PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));CHKERRQ(ierr);
110f268cba8SShri Abhyankar       array += nnz;
111f268cba8SShri Abhyankar       row++;
112f268cba8SShri Abhyankar     }
113f268cba8SShri Abhyankar 
114f268cba8SShri Abhyankar     a     = (*head)->more_space;
115f268cba8SShri Abhyankar     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
116f268cba8SShri Abhyankar     ierr  = PetscFree(*head);CHKERRQ(ierr);
117f268cba8SShri Abhyankar     *head = a;
118f268cba8SShri Abhyankar   }
119f268cba8SShri Abhyankar   bi[n] = bi[n-1] + nnzL;
120f268cba8SShri Abhyankar   bdiag[n] = bdiag[n-1]-1;
121f268cba8SShri Abhyankar   PetscFunctionReturn(0);
122f268cba8SShri Abhyankar }
123f268cba8SShri Abhyankar 
124783ef271SHong Zhang /*
125783ef271SHong Zhang   PetscFreeSpaceContiguous_Cholesky -
126783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array
127783ef271SHong Zhang   that enables an efficient matrix triangular solve.
128783ef271SHong Zhang 
129783ef271SHong Zhang    Input Parameters:
130783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization
131783ef271SHong Zhang .  space - an allocated int array with length nnz of factored matrix.
132783ef271SHong Zhang .  n - order of the matrix
133783ef271SHong Zhang .  ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix.
134783ef271SHong Zhang -  udiag - int array of length n.
135783ef271SHong Zhang 
136783ef271SHong Zhang    Output Parameter:
137783ef271SHong Zhang +  space - column indices are copied into this int array with contiguous layout of U, with diagonal located as the last entry in each row
138783ef271SHong Zhang -  udiag - indices of diagonal entries
139783ef271SHong Zhang 
140783ef271SHong Zhang    See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description.
141783ef271SHong Zhang */
142783ef271SHong Zhang 
143783ef271SHong Zhang #undef __FUNCT__
144783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_Cholesky"
145783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *ui,PetscInt *udiag)
146783ef271SHong Zhang {
147783ef271SHong Zhang   PetscFreeSpaceList a;
148783ef271SHong Zhang   PetscErrorCode     ierr;
149783ef271SHong Zhang   PetscInt           row,nnz,*uj,*array,total;
150783ef271SHong Zhang 
151783ef271SHong Zhang   PetscFunctionBegin;
152783ef271SHong Zhang   row   = 0;
153783ef271SHong Zhang   total = 0;
154783ef271SHong Zhang   while ((*head)!=NULL) {
155783ef271SHong Zhang     total += (*head)->local_used;
156783ef271SHong Zhang     array  = (*head)->array_head;
157783ef271SHong Zhang 
158783ef271SHong Zhang     while (ui[row+1] <= total && row < n){
159783ef271SHong Zhang       udiag[row] = ui[row+1] - 1;     /* points to the last entry of U(row,:) */
160783ef271SHong Zhang       nnz  = ui[row+1] - ui[row] - 1; /* exclude diagonal */
161783ef271SHong Zhang       uj   = space + ui[row];
162783ef271SHong Zhang       ierr = PetscMemcpy(uj,array+1,nnz*sizeof(PetscInt));CHKERRQ(ierr);
163783ef271SHong Zhang       uj[nnz] = array[0]; /* diagonal */
164783ef271SHong Zhang       array += nnz + 1;
165783ef271SHong Zhang       row++;
166783ef271SHong Zhang     }
167783ef271SHong Zhang 
168783ef271SHong Zhang     a     = (*head)->more_space;
169783ef271SHong Zhang     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
170783ef271SHong Zhang     ierr  = PetscFree(*head);CHKERRQ(ierr);
171783ef271SHong Zhang     *head = a;
172783ef271SHong Zhang   }
173783ef271SHong Zhang   PetscFunctionReturn(0);
174783ef271SHong Zhang }
175783ef271SHong Zhang 
17612b5cbf3SHong Zhang #undef __FUNCT__
177a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceDestroy"
178a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
1797a48dd6fSHong Zhang {
180a1a86e44SBarry Smith   PetscFreeSpaceList a;
1817a48dd6fSHong Zhang   PetscErrorCode     ierr;
1827a48dd6fSHong Zhang 
1837a48dd6fSHong Zhang   PetscFunctionBegin;
1847a48dd6fSHong Zhang   while ((head)!=NULL) {
1857a48dd6fSHong Zhang     a    = (head)->more_space;
1867a48dd6fSHong Zhang     ierr = PetscFree((head)->array_head);CHKERRQ(ierr);
1877a48dd6fSHong Zhang     ierr = PetscFree(head);CHKERRQ(ierr);
1887a48dd6fSHong Zhang     head = a;
1897a48dd6fSHong Zhang   }
1907a48dd6fSHong Zhang   PetscFunctionReturn(0);
1917a48dd6fSHong Zhang }
192