xref: /petsc/src/mat/utils/freespace.c (revision 783ef271488b82261007e3dc5c21163b50ca651c)
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;
3970f19b1fSKris Buschelman   while ((*head)!=NULL) {
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 /*
51*783ef271SHong Zhang   PetscFreeSpaceContiguous_LU -
52*783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array
53*783ef271SHong Zhang   that enables an efficient matrix triangular solve.
5430cb48eeSHong Zhang 
5530cb48eeSHong Zhang    Input Parameters:
56*783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ILU or LU factorization
57*783ef271SHong Zhang .  space - an allocated int array with length nnz of factored matrix.
5830cb48eeSHong Zhang .  n - order of the matrix
5930cb48eeSHong Zhang .  bi - row pointer of factored matrix with length 2n+2. First n+1 entries are set based on the traditional layout of L and U matrices
6030cb48eeSHong Zhang -  bdiag - int array holding the number of nonzeros in each row of L matrix, excluding diagonals.
6130cb48eeSHong Zhang 
6230cb48eeSHong Zhang    Output Parameter:
6330cb48eeSHong Zhang .  space - column indices are copied into this int array with contiguous layout of L and U
64*783ef271SHong Zhang 
6530cb48eeSHong Zhang    See MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct() for detailed description.
6630cb48eeSHong Zhang */
677a48dd6fSHong Zhang #undef __FUNCT__
68*783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_LU"
69*783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
7012b5cbf3SHong Zhang {
7112b5cbf3SHong Zhang   PetscFreeSpaceList a;
7212b5cbf3SHong Zhang   PetscErrorCode     ierr;
7312b5cbf3SHong Zhang   PetscInt           row,nnz,*bj,*array,total;
7412b5cbf3SHong Zhang   PetscInt           nnzL,nnzU;
7512b5cbf3SHong Zhang 
7612b5cbf3SHong Zhang   PetscFunctionBegin;
7712b5cbf3SHong Zhang   bi[2*n+1] = bi[n];
7830cb48eeSHong Zhang   row       = 0;
7912b5cbf3SHong Zhang   total     = 0;
8012b5cbf3SHong Zhang   nnzL  = bdiag[0];
8112b5cbf3SHong Zhang   while ((*head)!=NULL) {
8212b5cbf3SHong Zhang     total += (*head)->local_used;
8312b5cbf3SHong Zhang     array  = (*head)->array_head;
8412b5cbf3SHong Zhang 
8530cb48eeSHong Zhang     while (bi[row+1] <= total && row < n){
8630cb48eeSHong Zhang       /* copy array entries into bj for this row */
8730cb48eeSHong Zhang       nnz  = bi[row+1] - bi[row];
8830cb48eeSHong Zhang       /* set bi[row] for new datastruct */
8930cb48eeSHong Zhang       if (row == 0 ){
9030cb48eeSHong Zhang         bi[row] = 0;
9112b5cbf3SHong Zhang       } else {
9230cb48eeSHong Zhang         bi[row] = bi[row-1] + nnzL; /* nnzL of previous row */
9312b5cbf3SHong Zhang       }
9412b5cbf3SHong Zhang 
9512b5cbf3SHong Zhang       /* L part */
9630cb48eeSHong Zhang       nnzL = bdiag[row];
9730cb48eeSHong Zhang       bj   = space+bi[row];
9812b5cbf3SHong Zhang       ierr = PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));CHKERRQ(ierr);
9912b5cbf3SHong Zhang 
10012b5cbf3SHong Zhang       /* diagonal entry */
10130cb48eeSHong Zhang       bdiag[row]        = bi[2*n-row+1]-1;
10230cb48eeSHong Zhang       space[bdiag[row]] = row;
10312b5cbf3SHong Zhang 
10412b5cbf3SHong Zhang       /* U part */
10512b5cbf3SHong Zhang       nnzU        = nnz - nnzL;
10630cb48eeSHong Zhang       bi[2*n-row] = bi[2*n-row+1] - nnzU;
10712b5cbf3SHong Zhang       nnzU --;      /* exclude diagonal */
10830cb48eeSHong Zhang       bj   = space + bi[2*n-(row)];
10912b5cbf3SHong Zhang       ierr = PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));CHKERRQ(ierr);
11012b5cbf3SHong Zhang 
11112b5cbf3SHong Zhang       array += nnz;
11212b5cbf3SHong Zhang       row++;
11312b5cbf3SHong Zhang     }
11412b5cbf3SHong Zhang 
11512b5cbf3SHong Zhang     a     = (*head)->more_space;
11612b5cbf3SHong Zhang     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
11712b5cbf3SHong Zhang     ierr  = PetscFree(*head);CHKERRQ(ierr);
11812b5cbf3SHong Zhang     *head = a;
11912b5cbf3SHong Zhang   }
12012b5cbf3SHong Zhang   bi[n] = bi[n-1] + nnzL;
12112b5cbf3SHong Zhang   if (bi[n] != bi[n+1]) SETERRQ2(1,"bi[n] %d != bi[n+1] %d",bi[n],bi[n+1]);
12212b5cbf3SHong Zhang   PetscFunctionReturn(0);
12312b5cbf3SHong Zhang }
12412b5cbf3SHong Zhang 
125*783ef271SHong Zhang /*
126*783ef271SHong Zhang   PetscFreeSpaceContiguous_Cholesky -
127*783ef271SHong Zhang     Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array
128*783ef271SHong Zhang   that enables an efficient matrix triangular solve.
129*783ef271SHong Zhang 
130*783ef271SHong Zhang    Input Parameters:
131*783ef271SHong Zhang +  head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization
132*783ef271SHong Zhang .  space - an allocated int array with length nnz of factored matrix.
133*783ef271SHong Zhang .  n - order of the matrix
134*783ef271SHong Zhang .  ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix.
135*783ef271SHong Zhang -  udiag - int array of length n.
136*783ef271SHong Zhang 
137*783ef271SHong Zhang    Output Parameter:
138*783ef271SHong 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
139*783ef271SHong Zhang -  udiag - indices of diagonal entries
140*783ef271SHong Zhang 
141*783ef271SHong Zhang    See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description.
142*783ef271SHong Zhang */
143*783ef271SHong Zhang 
144*783ef271SHong Zhang #undef __FUNCT__
145*783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_Cholesky"
146*783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *ui,PetscInt *udiag)
147*783ef271SHong Zhang {
148*783ef271SHong Zhang   PetscFreeSpaceList a;
149*783ef271SHong Zhang   PetscErrorCode     ierr;
150*783ef271SHong Zhang   PetscInt           row,nnz,*uj,*array,total;
151*783ef271SHong Zhang 
152*783ef271SHong Zhang   PetscFunctionBegin;
153*783ef271SHong Zhang   row   = 0;
154*783ef271SHong Zhang   total = 0;
155*783ef271SHong Zhang   while ((*head)!=NULL) {
156*783ef271SHong Zhang     total += (*head)->local_used;
157*783ef271SHong Zhang     array  = (*head)->array_head;
158*783ef271SHong Zhang 
159*783ef271SHong Zhang     while (ui[row+1] <= total && row < n){
160*783ef271SHong Zhang       udiag[row] = ui[row+1] - 1;     /* points to the last entry of U(row,:) */
161*783ef271SHong Zhang       nnz  = ui[row+1] - ui[row] - 1; /* exclude diagonal */
162*783ef271SHong Zhang       uj   = space + ui[row];
163*783ef271SHong Zhang       ierr = PetscMemcpy(uj,array+1,nnz*sizeof(PetscInt));CHKERRQ(ierr);
164*783ef271SHong Zhang       uj[nnz] = array[0]; /* diagonal */
165*783ef271SHong Zhang       array += nnz + 1;
166*783ef271SHong Zhang       row++;
167*783ef271SHong Zhang     }
168*783ef271SHong Zhang 
169*783ef271SHong Zhang     a     = (*head)->more_space;
170*783ef271SHong Zhang     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
171*783ef271SHong Zhang     ierr  = PetscFree(*head);CHKERRQ(ierr);
172*783ef271SHong Zhang     *head = a;
173*783ef271SHong Zhang   }
174*783ef271SHong Zhang   PetscFunctionReturn(0);
175*783ef271SHong Zhang }
176*783ef271SHong Zhang 
17712b5cbf3SHong Zhang #undef __FUNCT__
178a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceDestroy"
179a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
1807a48dd6fSHong Zhang {
181a1a86e44SBarry Smith   PetscFreeSpaceList a;
1827a48dd6fSHong Zhang   PetscErrorCode     ierr;
1837a48dd6fSHong Zhang 
1847a48dd6fSHong Zhang   PetscFunctionBegin;
1857a48dd6fSHong Zhang   while ((head)!=NULL) {
1867a48dd6fSHong Zhang     a    = (head)->more_space;
1877a48dd6fSHong Zhang     ierr = PetscFree((head)->array_head);CHKERRQ(ierr);
1887a48dd6fSHong Zhang     ierr = PetscFree(head);CHKERRQ(ierr);
1897a48dd6fSHong Zhang     head = a;
1907a48dd6fSHong Zhang   }
1917a48dd6fSHong Zhang   PetscFunctionReturn(0);
1927a48dd6fSHong Zhang }
193