1be1d678aSKris Buschelman 283287d42SBarry Smith /* 383287d42SBarry Smith Factorization code for BAIJ format. 483287d42SBarry Smith */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 6c6db04a5SJed Brown #include <../src/mat/blockinvert.h> 783287d42SBarry Smith 8db4efbfdSBarry Smith #undef __FUNCT__ 94dd39f65SShri Abhyankar #define __FUNCT__ "MatSeqBAIJSetNumericFactorization" 10db4efbfdSBarry Smith /* 11db4efbfdSBarry Smith This is used to set the numeric factorization for both LU and ILU symbolic factorization 12db4efbfdSBarry Smith */ 13ace3abfcSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat fact,PetscBool natural) 14ae3d28f0SHong Zhang { 15ae3d28f0SHong Zhang PetscFunctionBegin; 166929473cSShri Abhyankar if(natural){ 176929473cSShri Abhyankar switch (fact->rmap->bs){ 18048b5e81SShri Abhyankar case 1: 19048b5e81SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 20048b5e81SShri Abhyankar break; 216929473cSShri Abhyankar case 2: 224dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 236929473cSShri Abhyankar break; 246929473cSShri Abhyankar case 3: 254dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 266929473cSShri Abhyankar break; 276929473cSShri Abhyankar case 4: 284dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 296929473cSShri Abhyankar break; 306929473cSShri Abhyankar case 5: 314dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 326929473cSShri Abhyankar break; 336929473cSShri Abhyankar case 6: 344dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 356929473cSShri Abhyankar break; 366929473cSShri Abhyankar case 7: 374dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 386929473cSShri Abhyankar break; 3929a97285SShri Abhyankar case 15: 4029a97285SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering; 4129a97285SShri Abhyankar break; 426929473cSShri Abhyankar default: 434dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 446929473cSShri Abhyankar break; 456929473cSShri Abhyankar } 466ba06ab7SHong Zhang } else { 47ae3d28f0SHong Zhang switch (fact->rmap->bs){ 48048b5e81SShri Abhyankar case 1: 49048b5e81SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 50048b5e81SShri Abhyankar break; 51ae3d28f0SHong Zhang case 2: 524dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 53ae3d28f0SHong Zhang break; 54ae3d28f0SHong Zhang case 3: 554dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 56ae3d28f0SHong Zhang break; 57ae3d28f0SHong Zhang case 4: 584dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 59ae3d28f0SHong Zhang break; 60ae3d28f0SHong Zhang case 5: 614dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 62ae3d28f0SHong Zhang break; 63ae3d28f0SHong Zhang case 6: 644dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 65ae3d28f0SHong Zhang break; 66ae3d28f0SHong Zhang case 7: 674dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7; 68ae3d28f0SHong Zhang break; 69ae3d28f0SHong Zhang default: 704dd39f65SShri Abhyankar fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 71ae3d28f0SHong Zhang break; 72ae3d28f0SHong Zhang } 736929473cSShri Abhyankar } 74ae3d28f0SHong Zhang PetscFunctionReturn(0); 75ae3d28f0SHong Zhang } 76ae3d28f0SHong Zhang 771008731dSBarry Smith #undef __FUNCT__ 781008731dSBarry Smith #define __FUNCT__ "MatSeqBAIJSetNumericFactorization_inplace" 79ace3abfcSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization_inplace(Mat inA,PetscBool natural) 80db4efbfdSBarry Smith { 81db4efbfdSBarry Smith PetscFunctionBegin; 82db4efbfdSBarry Smith if (natural) { 83db4efbfdSBarry Smith switch (inA->rmap->bs) { 84db4efbfdSBarry Smith case 1: 8506e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace; 86db4efbfdSBarry Smith break; 87db4efbfdSBarry Smith case 2: 8806e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace; 89db4efbfdSBarry Smith break; 90db4efbfdSBarry Smith case 3: 9106e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace; 92db4efbfdSBarry Smith break; 93db4efbfdSBarry Smith case 4: 94*ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE) 95db4efbfdSBarry Smith { 96ace3abfcSBarry Smith PetscBool sse_enabled_local; 97db4efbfdSBarry Smith PetscErrorCode ierr; 98db4efbfdSBarry Smith ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,PETSC_NULL);CHKERRQ(ierr); 99db4efbfdSBarry Smith if (sse_enabled_local) { 100db4efbfdSBarry Smith # if defined(PETSC_HAVE_SSE) 101db4efbfdSBarry Smith int i,*AJ=a->j,nz=a->nz,n=a->mbs; 102db4efbfdSBarry Smith if (n==(unsigned short)n) { 103db4efbfdSBarry Smith unsigned short *aj=(unsigned short *)AJ; 104db4efbfdSBarry Smith for (i=0;i<nz;i++) { 105db4efbfdSBarry Smith aj[i] = (unsigned short)AJ[i]; 106db4efbfdSBarry Smith } 107db4efbfdSBarry Smith inA->ops->setunfactored = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj; 108db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj; 109db4efbfdSBarry Smith ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr); 110db4efbfdSBarry Smith } else { 111db4efbfdSBarry Smith /* Scale the column indices for easier indexing in MatSolve. */ 112db4efbfdSBarry Smith /* for (i=0;i<nz;i++) { */ 113db4efbfdSBarry Smith /* AJ[i] = AJ[i]*4; */ 114db4efbfdSBarry Smith /* } */ 115db4efbfdSBarry Smith inA->ops->setunfactored = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE; 116db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE; 117db4efbfdSBarry Smith ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr); 118db4efbfdSBarry Smith } 119db4efbfdSBarry Smith # else 120db4efbfdSBarry Smith /* This should never be reached. If so, problem in PetscSSEIsEnabled. */ 121e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSE Hardware unavailable"); 122db4efbfdSBarry Smith # endif 123db4efbfdSBarry Smith } else { 12406e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace; 125db4efbfdSBarry Smith } 126db4efbfdSBarry Smith } 127db4efbfdSBarry Smith #else 12806e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace; 129db4efbfdSBarry Smith #endif 130db4efbfdSBarry Smith break; 131db4efbfdSBarry Smith case 5: 13206e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace; 133db4efbfdSBarry Smith break; 134db4efbfdSBarry Smith case 6: 13506e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace; 136db4efbfdSBarry Smith break; 137db4efbfdSBarry Smith case 7: 13806e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace; 139db4efbfdSBarry Smith break; 140db4efbfdSBarry Smith default: 14106e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace; 142db4efbfdSBarry Smith break; 143db4efbfdSBarry Smith } 144db4efbfdSBarry Smith } else { 145db4efbfdSBarry Smith switch (inA->rmap->bs) { 146db4efbfdSBarry Smith case 1: 14706e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace; 148db4efbfdSBarry Smith break; 149db4efbfdSBarry Smith case 2: 15006e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace; 151db4efbfdSBarry Smith break; 152db4efbfdSBarry Smith case 3: 15306e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace; 154db4efbfdSBarry Smith break; 155db4efbfdSBarry Smith case 4: 15606e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace; 157db4efbfdSBarry Smith break; 158db4efbfdSBarry Smith case 5: 15906e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace; 160db4efbfdSBarry Smith break; 161db4efbfdSBarry Smith case 6: 16206e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace; 163db4efbfdSBarry Smith break; 164db4efbfdSBarry Smith case 7: 16506e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace; 166db4efbfdSBarry Smith break; 167db4efbfdSBarry Smith default: 16806e38f1dSHong Zhang inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace; 169db4efbfdSBarry Smith break; 170db4efbfdSBarry Smith } 171db4efbfdSBarry Smith } 172db4efbfdSBarry Smith PetscFunctionReturn(0); 173db4efbfdSBarry Smith } 174db4efbfdSBarry Smith 17583287d42SBarry Smith /* 17683287d42SBarry Smith The symbolic factorization code is identical to that for AIJ format, 17783287d42SBarry Smith except for very small changes since this is now a SeqBAIJ datastructure. 17883287d42SBarry Smith NOT good code reuse. 17983287d42SBarry Smith */ 180c6db04a5SJed Brown #include <petscbt.h> 181c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 182c7272abeSHong Zhang 1834a2ae208SSatish Balay #undef __FUNCT__ 1844dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ" 1854dd39f65SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 18683287d42SBarry Smith { 18783287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 188c7272abeSHong Zhang PetscInt n=a->mbs,bs = A->rmap->bs,bs2=a->bs2; 189ace3abfcSBarry Smith PetscBool row_identity,col_identity,both_identity; 19083287d42SBarry Smith IS isicol; 1916849ba73SBarry Smith PetscErrorCode ierr; 192c7272abeSHong Zhang const PetscInt *r,*ic; 193c7272abeSHong Zhang PetscInt i,*ai=a->i,*aj=a->j; 194c7272abeSHong Zhang PetscInt *bi,*bj,*ajtmp; 195c7272abeSHong Zhang PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 196c7272abeSHong Zhang PetscReal f; 197c7272abeSHong Zhang PetscInt nlnk,*lnk,k,**bi_ptr; 198c7272abeSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 199c7272abeSHong Zhang PetscBT lnkbt; 20083287d42SBarry Smith 20183287d42SBarry Smith PetscFunctionBegin; 202e32f2f54SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square"); 2036ba06ab7SHong Zhang if (bs>1){ /* check shifttype */ 2046ba06ab7SHong Zhang if (info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) 2056ba06ab7SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix"); 2066ba06ab7SHong Zhang } 2076ba06ab7SHong Zhang 20883287d42SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 20983287d42SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 21083287d42SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 21183287d42SBarry Smith 212bc74c81fSJed Brown /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */ 213bc74c81fSJed Brown ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 214bc74c81fSJed Brown ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 215a6df321fSShri Abhyankar bi[0] = bdiag[0] = 0; 216b2b2dd24SShri Abhyankar 217b2b2dd24SShri Abhyankar /* linked list for storing column indices of the active row */ 218b2b2dd24SShri Abhyankar nlnk = n + 1; 219b2b2dd24SShri Abhyankar ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 220b2b2dd24SShri Abhyankar 221b2b2dd24SShri Abhyankar ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 222b2b2dd24SShri Abhyankar 223b2b2dd24SShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 224b2b2dd24SShri Abhyankar f = info->fill; 225b2b2dd24SShri Abhyankar ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 226b2b2dd24SShri Abhyankar current_space = free_space; 227b2b2dd24SShri Abhyankar 228b2b2dd24SShri Abhyankar for (i=0; i<n; i++) { 229b2b2dd24SShri Abhyankar /* copy previous fill into linked list */ 230b2b2dd24SShri Abhyankar nzi = 0; 231b2b2dd24SShri Abhyankar nnz = ai[r[i]+1] - ai[r[i]]; 232e32f2f54SBarry Smith if (!nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 233b2b2dd24SShri Abhyankar ajtmp = aj + ai[r[i]]; 234b2b2dd24SShri Abhyankar ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 235b2b2dd24SShri Abhyankar nzi += nlnk; 236b2b2dd24SShri Abhyankar 237b2b2dd24SShri Abhyankar /* add pivot rows into linked list */ 238b2b2dd24SShri Abhyankar row = lnk[n]; 239b2b2dd24SShri Abhyankar while (row < i) { 240b2b2dd24SShri Abhyankar nzbd = bdiag[row] + 1; /* num of entries in the row with column index <= row */ 241b2b2dd24SShri Abhyankar ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 242b2b2dd24SShri Abhyankar ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 243b2b2dd24SShri Abhyankar nzi += nlnk; 244b2b2dd24SShri Abhyankar row = lnk[row]; 245b2b2dd24SShri Abhyankar } 246b2b2dd24SShri Abhyankar bi[i+1] = bi[i] + nzi; 247b2b2dd24SShri Abhyankar im[i] = nzi; 248b2b2dd24SShri Abhyankar 249b2b2dd24SShri Abhyankar /* mark bdiag */ 250b2b2dd24SShri Abhyankar nzbd = 0; 251b2b2dd24SShri Abhyankar nnz = nzi; 252b2b2dd24SShri Abhyankar k = lnk[n]; 253b2b2dd24SShri Abhyankar while (nnz-- && k < i){ 254b2b2dd24SShri Abhyankar nzbd++; 255b2b2dd24SShri Abhyankar k = lnk[k]; 256b2b2dd24SShri Abhyankar } 2572ce24eb6SHong Zhang bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */ 258b2b2dd24SShri Abhyankar 259b2b2dd24SShri Abhyankar /* if free space is not available, make more free space */ 260b2b2dd24SShri Abhyankar if (current_space->local_remaining<nzi) { 261b2b2dd24SShri Abhyankar nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */ 262b2b2dd24SShri Abhyankar ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 263b2b2dd24SShri Abhyankar reallocs++; 264b2b2dd24SShri Abhyankar } 265b2b2dd24SShri Abhyankar 266b2b2dd24SShri Abhyankar /* copy data into free space, then initialize lnk */ 267b2b2dd24SShri Abhyankar ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 268b2b2dd24SShri Abhyankar bi_ptr[i] = current_space->array; 269b2b2dd24SShri Abhyankar current_space->array += nzi; 270b2b2dd24SShri Abhyankar current_space->local_used += nzi; 271b2b2dd24SShri Abhyankar current_space->local_remaining -= nzi; 272b2b2dd24SShri Abhyankar } 273b2b2dd24SShri Abhyankar 274b2b2dd24SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 275b2b2dd24SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 276b2b2dd24SShri Abhyankar 2779263d837SHong Zhang /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 278b2b2dd24SShri Abhyankar ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 2792ce24eb6SHong Zhang ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 280b2b2dd24SShri Abhyankar ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 281b2b2dd24SShri Abhyankar ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 282b2b2dd24SShri Abhyankar 283b2b2dd24SShri Abhyankar /* put together the new matrix */ 284b2b2dd24SShri Abhyankar ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 285b2b2dd24SShri Abhyankar ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 286b2b2dd24SShri Abhyankar b = (Mat_SeqBAIJ*)(B)->data; 287b2b2dd24SShri Abhyankar b->free_a = PETSC_TRUE; 288b2b2dd24SShri Abhyankar b->free_ij = PETSC_TRUE; 289b2b2dd24SShri Abhyankar b->singlemalloc = PETSC_FALSE; 290b2b2dd24SShri Abhyankar ierr = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 291b2b2dd24SShri Abhyankar b->j = bj; 292b2b2dd24SShri Abhyankar b->i = bi; 293b2b2dd24SShri Abhyankar b->diag = bdiag; 294b2b2dd24SShri Abhyankar b->free_diag = PETSC_TRUE; 295b2b2dd24SShri Abhyankar b->ilen = 0; 296b2b2dd24SShri Abhyankar b->imax = 0; 297b2b2dd24SShri Abhyankar b->row = isrow; 298b2b2dd24SShri Abhyankar b->col = iscol; 299b2b2dd24SShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 300b2b2dd24SShri Abhyankar ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 301b2b2dd24SShri Abhyankar ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 302b2b2dd24SShri Abhyankar b->icol = isicol; 303b2b2dd24SShri Abhyankar ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 304b2b2dd24SShri Abhyankar ierr = PetscLogObjectMemory(B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr); 305b2b2dd24SShri Abhyankar 306b2b2dd24SShri Abhyankar b->maxnz = b->nz = bdiag[0]+1; 307d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 308b2b2dd24SShri Abhyankar B->info.factor_mallocs = reallocs; 309b2b2dd24SShri Abhyankar B->info.fill_ratio_given = f; 310b2b2dd24SShri Abhyankar 311b2b2dd24SShri Abhyankar if (ai[n] != 0) { 312b2b2dd24SShri Abhyankar B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 313b2b2dd24SShri Abhyankar } else { 314b2b2dd24SShri Abhyankar B->info.fill_ratio_needed = 0.0; 315b2b2dd24SShri Abhyankar } 3169263d837SHong Zhang #if defined(PETSC_USE_INFO) 3179263d837SHong Zhang if (ai[n] != 0) { 3189263d837SHong Zhang PetscReal af = B->info.fill_ratio_needed; 3199263d837SHong Zhang ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 3209263d837SHong Zhang ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 3219263d837SHong Zhang ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 3229263d837SHong Zhang ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 3239263d837SHong Zhang } else { 3249263d837SHong Zhang ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 3259263d837SHong Zhang } 3269263d837SHong Zhang #endif 327b2b2dd24SShri Abhyankar 328b2b2dd24SShri Abhyankar ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 329b2b2dd24SShri Abhyankar ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 330ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 3314dd39f65SShri Abhyankar ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr); 332b2b2dd24SShri Abhyankar PetscFunctionReturn(0); 333b2b2dd24SShri Abhyankar } 334b2b2dd24SShri Abhyankar 335b2b2dd24SShri Abhyankar #undef __FUNCT__ 33606e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_inplace" 33706e38f1dSHong Zhang PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 338faca2338SShri Abhyankar { 339faca2338SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 340faca2338SShri Abhyankar PetscInt n=a->mbs,bs = A->rmap->bs,bs2=a->bs2; 341ace3abfcSBarry Smith PetscBool row_identity,col_identity,both_identity; 342faca2338SShri Abhyankar IS isicol; 343faca2338SShri Abhyankar PetscErrorCode ierr; 344faca2338SShri Abhyankar const PetscInt *r,*ic; 345faca2338SShri Abhyankar PetscInt i,*ai=a->i,*aj=a->j; 346faca2338SShri Abhyankar PetscInt *bi,*bj,*ajtmp; 347faca2338SShri Abhyankar PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 348faca2338SShri Abhyankar PetscReal f; 349faca2338SShri Abhyankar PetscInt nlnk,*lnk,k,**bi_ptr; 350faca2338SShri Abhyankar PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 351faca2338SShri Abhyankar PetscBT lnkbt; 352faca2338SShri Abhyankar 353faca2338SShri Abhyankar PetscFunctionBegin; 354e32f2f54SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square"); 355faca2338SShri Abhyankar ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 356faca2338SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 357faca2338SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 358faca2338SShri Abhyankar 359bc74c81fSJed Brown /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */ 360bc74c81fSJed Brown ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 361bc74c81fSJed Brown ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 362bc74c81fSJed Brown 363a6df321fSShri Abhyankar bi[0] = bdiag[0] = 0; 364c7272abeSHong Zhang 365c7272abeSHong Zhang /* linked list for storing column indices of the active row */ 366c7272abeSHong Zhang nlnk = n + 1; 367c7272abeSHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 368c7272abeSHong Zhang 369c7272abeSHong Zhang ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 370c7272abeSHong Zhang 371c7272abeSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 372c7272abeSHong Zhang f = info->fill; 373c7272abeSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 374c7272abeSHong Zhang current_space = free_space; 37583287d42SBarry Smith 37683287d42SBarry Smith for (i=0; i<n; i++) { 377c7272abeSHong Zhang /* copy previous fill into linked list */ 37883287d42SBarry Smith nzi = 0; 379c7272abeSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 380e32f2f54SBarry Smith if (!nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 381c7272abeSHong Zhang ajtmp = aj + ai[r[i]]; 382c7272abeSHong Zhang ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 383c7272abeSHong Zhang nzi += nlnk; 384c7272abeSHong Zhang 385c7272abeSHong Zhang /* add pivot rows into linked list */ 386c7272abeSHong Zhang row = lnk[n]; 387c7272abeSHong Zhang while (row < i) { 388c7272abeSHong Zhang nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 389c7272abeSHong Zhang ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 390c7272abeSHong Zhang ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 391c7272abeSHong Zhang nzi += nlnk; 392c7272abeSHong Zhang row = lnk[row]; 39383287d42SBarry Smith } 394c7272abeSHong Zhang bi[i+1] = bi[i] + nzi; 395c7272abeSHong Zhang im[i] = nzi; 396c7272abeSHong Zhang 397c7272abeSHong Zhang /* mark bdiag */ 398c7272abeSHong Zhang nzbd = 0; 399c7272abeSHong Zhang nnz = nzi; 400c7272abeSHong Zhang k = lnk[n]; 401c7272abeSHong Zhang while (nnz-- && k < i){ 402c7272abeSHong Zhang nzbd++; 403c7272abeSHong Zhang k = lnk[k]; 404c7272abeSHong Zhang } 405c7272abeSHong Zhang bdiag[i] = bi[i] + nzbd; 406c7272abeSHong Zhang 407c7272abeSHong Zhang /* if free space is not available, make more free space */ 408c7272abeSHong Zhang if (current_space->local_remaining<nzi) { 409c7272abeSHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 410c7272abeSHong Zhang ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 411c7272abeSHong Zhang reallocs++; 41283287d42SBarry Smith } 41383287d42SBarry Smith 414c7272abeSHong Zhang /* copy data into free space, then initialize lnk */ 415c7272abeSHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 416c7272abeSHong Zhang bi_ptr[i] = current_space->array; 417c7272abeSHong Zhang current_space->array += nzi; 418c7272abeSHong Zhang current_space->local_used += nzi; 419c7272abeSHong Zhang current_space->local_remaining -= nzi; 420c7272abeSHong Zhang } 4216cf91177SBarry Smith #if defined(PETSC_USE_INFO) 42283287d42SBarry Smith if (ai[n] != 0) { 423c7272abeSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 424ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 425ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 426ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 427ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 42883287d42SBarry Smith } else { 429c7272abeSHong Zhang ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 43083287d42SBarry Smith } 43163ba0a88SBarry Smith #endif 43283287d42SBarry Smith 43383287d42SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 43483287d42SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 43583287d42SBarry Smith 436c7272abeSHong Zhang /* destroy list of free space and other temporary array(s) */ 437c7272abeSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 438c7272abeSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 439c7272abeSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 440c7272abeSHong Zhang ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 44183287d42SBarry Smith 44283287d42SBarry Smith /* put together the new matrix */ 443719d5645SBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 444719d5645SBarry Smith ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 445719d5645SBarry Smith b = (Mat_SeqBAIJ*)(B)->data; 446e6b907acSBarry Smith b->free_a = PETSC_TRUE; 447e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 448c7272abeSHong Zhang b->singlemalloc = PETSC_FALSE; 449c7272abeSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 450c7272abeSHong Zhang b->j = bj; 451c7272abeSHong Zhang b->i = bi; 452c7272abeSHong Zhang b->diag = bdiag; 4537f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 45483287d42SBarry Smith b->ilen = 0; 45583287d42SBarry Smith b->imax = 0; 45683287d42SBarry Smith b->row = isrow; 45783287d42SBarry Smith b->col = iscol; 458d3d32019SBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 45983287d42SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 46083287d42SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 46183287d42SBarry Smith b->icol = isicol; 46287828ca2SBarry Smith ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 463c7272abeSHong Zhang ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr); 46483287d42SBarry Smith 465c7272abeSHong Zhang b->maxnz = b->nz = bi[n] ; 466d5f3da31SBarry Smith (B)->factortype = MAT_FACTOR_LU; 467719d5645SBarry Smith (B)->info.factor_mallocs = reallocs; 468719d5645SBarry Smith (B)->info.fill_ratio_given = f; 469c7272abeSHong Zhang 47083287d42SBarry Smith if (ai[n] != 0) { 471c7272abeSHong Zhang (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 47283287d42SBarry Smith } else { 473719d5645SBarry Smith (B)->info.fill_ratio_needed = 0.0; 47483287d42SBarry Smith } 475c7272abeSHong Zhang 476db4efbfdSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 477db4efbfdSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 478ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 4798b1456e3SHong Zhang ierr = MatSeqBAIJSetNumericFactorization_inplace(B,both_identity);CHKERRQ(ierr); 48083287d42SBarry Smith PetscFunctionReturn(0); 48183287d42SBarry Smith } 482db4efbfdSBarry Smith 483