1be1d678aSKris Buschelman 2c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 370f19b1fSKris Buschelman 470f19b1fSKris Buschelman #undef __FUNCT__ 5a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceGet" 6a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceGet(PetscInt n,PetscFreeSpaceList *list) 72e111b49SBarry Smith { 8a1a86e44SBarry Smith PetscFreeSpaceList a; 9dfbe8321SBarry Smith PetscErrorCode ierr; 1070f19b1fSKris Buschelman 1170f19b1fSKris Buschelman PetscFunctionBegin; 12a1a86e44SBarry Smith ierr = PetscMalloc(sizeof(struct _Space),&a);CHKERRQ(ierr); 132e111b49SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&(a->array_head));CHKERRQ(ierr); 1470f19b1fSKris Buschelman a->array = a->array_head; 152e111b49SBarry Smith a->local_remaining = n; 1670f19b1fSKris Buschelman a->local_used = 0; 1770f19b1fSKris Buschelman a->total_array_size = 0; 1870f19b1fSKris Buschelman a->more_space = NULL; 1970f19b1fSKris Buschelman 2070f19b1fSKris Buschelman if (*list) { 2170f19b1fSKris Buschelman (*list)->more_space = a; 2270f19b1fSKris Buschelman a->total_array_size = (*list)->total_array_size; 2370f19b1fSKris Buschelman } 2470f19b1fSKris Buschelman 252e111b49SBarry Smith a->total_array_size += n; 2670f19b1fSKris Buschelman *list = a; 2770f19b1fSKris Buschelman PetscFunctionReturn(0); 2870f19b1fSKris Buschelman } 2970f19b1fSKris Buschelman 3070f19b1fSKris Buschelman #undef __FUNCT__ 31a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceContiguous" 32a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head,PetscInt *space) 3338baddfdSBarry Smith { 34a1a86e44SBarry Smith PetscFreeSpaceList a; 35dfbe8321SBarry Smith PetscErrorCode ierr; 3670f19b1fSKris Buschelman 3770f19b1fSKris Buschelman PetscFunctionBegin; 38c05d87d6SBarry Smith while ((*head)) { 3970f19b1fSKris Buschelman a = (*head)->more_space; 402e111b49SBarry Smith ierr = PetscMemcpy(space,(*head)->array_head,((*head)->local_used)*sizeof(PetscInt));CHKERRQ(ierr); 4170f19b1fSKris Buschelman space += (*head)->local_used; 4270f19b1fSKris Buschelman ierr = PetscFree((*head)->array_head);CHKERRQ(ierr); 4370f19b1fSKris Buschelman ierr = PetscFree(*head);CHKERRQ(ierr); 4470f19b1fSKris Buschelman *head = a; 4570f19b1fSKris Buschelman } 4670f19b1fSKris Buschelman PetscFunctionReturn(0); 4770f19b1fSKris Buschelman } 487a48dd6fSHong Zhang 4930cb48eeSHong Zhang /* 50783ef271SHong Zhang PetscFreeSpaceContiguous_LU - 51783ef271SHong Zhang Copy a linket list obtained from matrix symbolic ILU or LU factorization into a contiguous array 52783ef271SHong Zhang that enables an efficient matrix triangular solve. 5330cb48eeSHong Zhang 5430cb48eeSHong Zhang Input Parameters: 55783ef271SHong Zhang + head - linked list of column indices obtained from matrix symbolic ILU or LU factorization 56783ef271SHong Zhang . space - an allocated int array with length nnz of factored matrix. 5730cb48eeSHong Zhang . n - order of the matrix 582ce24eb6SHong Zhang . bi - row pointer of factored matrix L with length n+1. 592ce24eb6SHong 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. 6030cb48eeSHong Zhang 6130cb48eeSHong Zhang Output Parameter: 6230cb48eeSHong Zhang . space - column indices are copied into this int array with contiguous layout of L and U 63783ef271SHong Zhang 642ce24eb6SHong Zhang See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed data structure of L and U 6530cb48eeSHong Zhang */ 667a48dd6fSHong Zhang #undef __FUNCT__ 67783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_LU" 68783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_LU(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag) 6912b5cbf3SHong Zhang { 7012b5cbf3SHong Zhang PetscFreeSpaceList a; 7112b5cbf3SHong Zhang PetscErrorCode ierr; 72f268cba8SShri Abhyankar PetscInt row,nnz,*bj,*array,total,bi_temp; 73f268cba8SShri Abhyankar PetscInt nnzL,nnzU; 74f268cba8SShri Abhyankar 75f268cba8SShri Abhyankar PetscFunctionBegin; 76f268cba8SShri Abhyankar bi_temp = bi[n]; 77f268cba8SShri Abhyankar row = 0; 78f268cba8SShri Abhyankar total = 0; 79f268cba8SShri Abhyankar nnzL = bdiag[0]; 80f268cba8SShri Abhyankar while ((*head)!=NULL) { 81f268cba8SShri Abhyankar total += (*head)->local_used; 82f268cba8SShri Abhyankar array = (*head)->array_head; 83f268cba8SShri Abhyankar 840e97c5f4SHong Zhang while (row < n) { 850e97c5f4SHong Zhang if (bi[row+1] > total) break; 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 } 119*43c97eaeSBarry Smith if (n) { 120f268cba8SShri Abhyankar bi[n] = bi[n-1] + nnzL; 121f268cba8SShri Abhyankar bdiag[n] = bdiag[n-1]-1; 122*43c97eaeSBarry Smith } 123f268cba8SShri Abhyankar PetscFunctionReturn(0); 124f268cba8SShri Abhyankar } 125f268cba8SShri Abhyankar 126783ef271SHong Zhang /* 127783ef271SHong Zhang PetscFreeSpaceContiguous_Cholesky - 128783ef271SHong Zhang Copy a linket list obtained from matrix symbolic ICC or Cholesky factorization into a contiguous array 129783ef271SHong Zhang that enables an efficient matrix triangular solve. 130783ef271SHong Zhang 131783ef271SHong Zhang Input Parameters: 132783ef271SHong Zhang + head - linked list of column indices obtained from matrix symbolic ICC or Cholesky factorization 133783ef271SHong Zhang . space - an allocated int array with length nnz of factored matrix. 134783ef271SHong Zhang . n - order of the matrix 135783ef271SHong Zhang . ui - row pointer of factored matrix with length n+1. All entries are set based on the traditional layout U matrix. 136783ef271SHong Zhang - udiag - int array of length n. 137783ef271SHong Zhang 138783ef271SHong Zhang Output Parameter: 139783ef271SHong 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 140783ef271SHong Zhang - udiag - indices of diagonal entries 141783ef271SHong Zhang 142783ef271SHong Zhang See MatICCFactorSymbolic_SeqAIJ_newdatastruct() for detailed description. 143783ef271SHong Zhang */ 144783ef271SHong Zhang 145783ef271SHong Zhang #undef __FUNCT__ 146783ef271SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_Cholesky" 147783ef271SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_Cholesky(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *ui,PetscInt *udiag) 148783ef271SHong Zhang { 149783ef271SHong Zhang PetscFreeSpaceList a; 150783ef271SHong Zhang PetscErrorCode ierr; 151783ef271SHong Zhang PetscInt row,nnz,*uj,*array,total; 152783ef271SHong Zhang 153783ef271SHong Zhang PetscFunctionBegin; 154783ef271SHong Zhang row = 0; 155783ef271SHong Zhang total = 0; 1560e97c5f4SHong Zhang while (*head) { 157783ef271SHong Zhang total += (*head)->local_used; 158783ef271SHong Zhang array = (*head)->array_head; 159783ef271SHong Zhang 1600e97c5f4SHong Zhang while (row < n){ 1610e97c5f4SHong Zhang if (ui[row+1] > total) break; 162783ef271SHong Zhang udiag[row] = ui[row+1] - 1; /* points to the last entry of U(row,:) */ 163783ef271SHong Zhang nnz = ui[row+1] - ui[row] - 1; /* exclude diagonal */ 164783ef271SHong Zhang uj = space + ui[row]; 165783ef271SHong Zhang ierr = PetscMemcpy(uj,array+1,nnz*sizeof(PetscInt));CHKERRQ(ierr); 166783ef271SHong Zhang uj[nnz] = array[0]; /* diagonal */ 167783ef271SHong Zhang array += nnz + 1; 168783ef271SHong Zhang row++; 169783ef271SHong Zhang } 170783ef271SHong Zhang 171783ef271SHong Zhang a = (*head)->more_space; 172783ef271SHong Zhang ierr = PetscFree((*head)->array_head);CHKERRQ(ierr); 173783ef271SHong Zhang ierr = PetscFree(*head);CHKERRQ(ierr); 174783ef271SHong Zhang *head = a; 175783ef271SHong Zhang } 176783ef271SHong Zhang PetscFunctionReturn(0); 177783ef271SHong Zhang } 178783ef271SHong Zhang 17912b5cbf3SHong Zhang #undef __FUNCT__ 180a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceDestroy" 181a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head) 1827a48dd6fSHong Zhang { 183a1a86e44SBarry Smith PetscFreeSpaceList a; 1847a48dd6fSHong Zhang PetscErrorCode ierr; 1857a48dd6fSHong Zhang 1867a48dd6fSHong Zhang PetscFunctionBegin; 1877a48dd6fSHong Zhang while ((head)!=NULL) { 1887a48dd6fSHong Zhang a = (head)->more_space; 1897a48dd6fSHong Zhang ierr = PetscFree((head)->array_head);CHKERRQ(ierr); 1907a48dd6fSHong Zhang ierr = PetscFree(head);CHKERRQ(ierr); 1917a48dd6fSHong Zhang head = a; 1927a48dd6fSHong Zhang } 1937a48dd6fSHong Zhang PetscFunctionReturn(0); 1947a48dd6fSHong Zhang } 195