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 /* 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 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 64783ef271SHong Zhang 6530cb48eeSHong Zhang See MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct() for detailed description. 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; 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*f268cba8SShri Abhyankar #undef __FUNCT__ 126*f268cba8SShri Abhyankar #define __FUNCT__ "PetscFreeSpaceContiguous_LU_v2" 127*f268cba8SShri Abhyankar PetscErrorCode PetscFreeSpaceContiguous_LU_v2(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag) 128*f268cba8SShri Abhyankar { 129*f268cba8SShri Abhyankar PetscFreeSpaceList a; 130*f268cba8SShri Abhyankar PetscErrorCode ierr; 131*f268cba8SShri Abhyankar PetscInt row,nnz,*bj,*array,total,bi_temp; 132*f268cba8SShri Abhyankar PetscInt nnzL,nnzU; 133*f268cba8SShri Abhyankar 134*f268cba8SShri Abhyankar PetscFunctionBegin; 135*f268cba8SShri Abhyankar /* bi[2*n+1] = bi[n]; */ 136*f268cba8SShri Abhyankar bi_temp = bi[n]; 137*f268cba8SShri Abhyankar row = 0; 138*f268cba8SShri Abhyankar total = 0; 139*f268cba8SShri Abhyankar nnzL = bdiag[0]; 140*f268cba8SShri Abhyankar while ((*head)!=NULL) { 141*f268cba8SShri Abhyankar total += (*head)->local_used; 142*f268cba8SShri Abhyankar array = (*head)->array_head; 143*f268cba8SShri Abhyankar 144*f268cba8SShri Abhyankar while (bi[row+1] <= total && row < n){ 145*f268cba8SShri Abhyankar /* copy array entries into bj for this row */ 146*f268cba8SShri Abhyankar nnz = bi[row+1] - bi[row]; 147*f268cba8SShri Abhyankar /* set bi[row] for new datastruct */ 148*f268cba8SShri Abhyankar if (row == 0 ){ 149*f268cba8SShri Abhyankar bi[row] = 0; 150*f268cba8SShri Abhyankar } else { 151*f268cba8SShri Abhyankar bi[row] = bi[row-1] + nnzL; /* nnzL of previous row */ 152*f268cba8SShri Abhyankar } 153*f268cba8SShri Abhyankar 154*f268cba8SShri Abhyankar /* L part */ 155*f268cba8SShri Abhyankar nnzL = bdiag[row]; 156*f268cba8SShri Abhyankar bj = space+bi[row]; 157*f268cba8SShri Abhyankar ierr = PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));CHKERRQ(ierr); 158*f268cba8SShri Abhyankar 159*f268cba8SShri Abhyankar /* diagonal entry */ 160*f268cba8SShri Abhyankar /* bdiag[row] = bi[2*n-row+1]-1; */ 161*f268cba8SShri Abhyankar bdiag[row] = bi_temp - 1; 162*f268cba8SShri Abhyankar space[bdiag[row]] = row; 163*f268cba8SShri Abhyankar 164*f268cba8SShri Abhyankar /* U part */ 165*f268cba8SShri Abhyankar nnzU = nnz - nnzL; 166*f268cba8SShri Abhyankar /* bi[2*n-row] = bi[2*n-row+1] - nnzU; */ 167*f268cba8SShri Abhyankar bi_temp = bi_temp - nnzU; 168*f268cba8SShri Abhyankar nnzU --; /* exclude diagonal */ 169*f268cba8SShri Abhyankar /* bj = space + bi[2*n-(row)]; */ 170*f268cba8SShri Abhyankar bj = space + bi_temp; 171*f268cba8SShri Abhyankar ierr = PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));CHKERRQ(ierr); 172*f268cba8SShri Abhyankar 173*f268cba8SShri Abhyankar array += nnz; 174*f268cba8SShri Abhyankar row++; 175*f268cba8SShri Abhyankar } 176*f268cba8SShri Abhyankar 177*f268cba8SShri Abhyankar a = (*head)->more_space; 178*f268cba8SShri Abhyankar ierr = PetscFree((*head)->array_head);CHKERRQ(ierr); 179*f268cba8SShri Abhyankar ierr = PetscFree(*head);CHKERRQ(ierr); 180*f268cba8SShri Abhyankar *head = a; 181*f268cba8SShri Abhyankar } 182*f268cba8SShri Abhyankar bi[n] = bi[n-1] + nnzL; 183*f268cba8SShri Abhyankar bdiag[n] = bdiag[n-1]-1; 184*f268cba8SShri Abhyankar /* if (bi[n] != bi[n+1]) SETERRQ2(1,"bi[n] %d != bi[n+1] %d",bi[n],bi[n+1]); */ 185*f268cba8SShri Abhyankar if(bi[n] != bdiag[n-1]) SETERRQ2(1,"bi[n] %d != bdiag[n-1] %d",bi[n],bdiag[n-1]); 186*f268cba8SShri Abhyankar PetscFunctionReturn(0); 187*f268cba8SShri Abhyankar } 188*f268cba8SShri Abhyankar 189*f268cba8SShri Abhyankar 190*f268cba8SShri Abhyankar 191783ef271SHong Zhang /* 192783ef271SHong Zhang PetscFreeSpaceContiguous_Cholesky - 193783ef271SHong Zhang Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array 194783ef271SHong Zhang that enables an efficient matrix triangular solve. 195783ef271SHong Zhang 196783ef271SHong Zhang Input Parameters: 197783ef271SHong Zhang + head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization 198783ef271SHong Zhang . space - an allocated int array with length nnz of factored matrix. 199783ef271SHong Zhang . n - order of the matrix 200783ef271SHong Zhang . ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix. 201783ef271SHong Zhang - udiag - int array of length n. 202783ef271SHong Zhang 203783ef271SHong Zhang Output Parameter: 204783ef271SHong 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 205783ef271SHong Zhang - udiag - indices of diagonal entries 206783ef271SHong Zhang 207783ef271SHong Zhang See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description. 208783ef271SHong Zhang */ 209783ef271SHong Zhang 210783ef271SHong Zhang #undef __FUNCT__ 211783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_Cholesky" 212783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *ui,PetscInt *udiag) 213783ef271SHong Zhang { 214783ef271SHong Zhang PetscFreeSpaceList a; 215783ef271SHong Zhang PetscErrorCode ierr; 216783ef271SHong Zhang PetscInt row,nnz,*uj,*array,total; 217783ef271SHong Zhang 218783ef271SHong Zhang PetscFunctionBegin; 219783ef271SHong Zhang row = 0; 220783ef271SHong Zhang total = 0; 221783ef271SHong Zhang while ((*head)!=NULL) { 222783ef271SHong Zhang total += (*head)->local_used; 223783ef271SHong Zhang array = (*head)->array_head; 224783ef271SHong Zhang 225783ef271SHong Zhang while (ui[row+1] <= total && row < n){ 226783ef271SHong Zhang udiag[row] = ui[row+1] - 1; /* points to the last entry of U(row,:) */ 227783ef271SHong Zhang nnz = ui[row+1] - ui[row] - 1; /* exclude diagonal */ 228783ef271SHong Zhang uj = space + ui[row]; 229783ef271SHong Zhang ierr = PetscMemcpy(uj,array+1,nnz*sizeof(PetscInt));CHKERRQ(ierr); 230783ef271SHong Zhang uj[nnz] = array[0]; /* diagonal */ 231783ef271SHong Zhang array += nnz + 1; 232783ef271SHong Zhang row++; 233783ef271SHong Zhang } 234783ef271SHong Zhang 235783ef271SHong Zhang a = (*head)->more_space; 236783ef271SHong Zhang ierr = PetscFree((*head)->array_head);CHKERRQ(ierr); 237783ef271SHong Zhang ierr = PetscFree(*head);CHKERRQ(ierr); 238783ef271SHong Zhang *head = a; 239783ef271SHong Zhang } 240783ef271SHong Zhang PetscFunctionReturn(0); 241783ef271SHong Zhang } 242783ef271SHong Zhang 24312b5cbf3SHong Zhang #undef __FUNCT__ 244a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceDestroy" 245a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head) 2467a48dd6fSHong Zhang { 247a1a86e44SBarry Smith PetscFreeSpaceList a; 2487a48dd6fSHong Zhang PetscErrorCode ierr; 2497a48dd6fSHong Zhang 2507a48dd6fSHong Zhang PetscFunctionBegin; 2517a48dd6fSHong Zhang while ((head)!=NULL) { 2527a48dd6fSHong Zhang a = (head)->more_space; 2537a48dd6fSHong Zhang ierr = PetscFree((head)->array_head);CHKERRQ(ierr); 2547a48dd6fSHong Zhang ierr = PetscFree(head);CHKERRQ(ierr); 2557a48dd6fSHong Zhang head = a; 2567a48dd6fSHong Zhang } 2577a48dd6fSHong Zhang PetscFunctionReturn(0); 2587a48dd6fSHong Zhang } 259