1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 383287d42SBarry Smith /* 483287d42SBarry Smith Factorization code for BAIJ format. 583287d42SBarry Smith */ 67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" 7c60f0209SBarry Smith #include "../src/mat/blockinvert.h" 883287d42SBarry Smith 9db4efbfdSBarry Smith #undef __FUNCT__ 10db4efbfdSBarry Smith #define __FUNCT__ "MatSeqBAIJSetNumericFactorization" 11db4efbfdSBarry Smith /* 12db4efbfdSBarry Smith This is used to set the numeric factorization for both LU and ILU symbolic factorization 13db4efbfdSBarry Smith */ 14db4efbfdSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat inA,PetscTruth natural) 15db4efbfdSBarry Smith { 16db4efbfdSBarry Smith PetscFunctionBegin; 17db4efbfdSBarry Smith if (natural) { 18db4efbfdSBarry Smith switch (inA->rmap->bs) { 19db4efbfdSBarry Smith case 1: 20db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 21db4efbfdSBarry Smith break; 22db4efbfdSBarry Smith case 2: 23db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 24db4efbfdSBarry Smith break; 25db4efbfdSBarry Smith case 3: 26db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 27db4efbfdSBarry Smith break; 28db4efbfdSBarry Smith case 4: 29db4efbfdSBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 30db4efbfdSBarry Smith { 31db4efbfdSBarry Smith PetscTruth sse_enabled_local; 32db4efbfdSBarry Smith PetscErrorCode ierr; 33db4efbfdSBarry Smith ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,PETSC_NULL);CHKERRQ(ierr); 34db4efbfdSBarry Smith if (sse_enabled_local) { 35db4efbfdSBarry Smith # if defined(PETSC_HAVE_SSE) 36db4efbfdSBarry Smith int i,*AJ=a->j,nz=a->nz,n=a->mbs; 37db4efbfdSBarry Smith if (n==(unsigned short)n) { 38db4efbfdSBarry Smith unsigned short *aj=(unsigned short *)AJ; 39db4efbfdSBarry Smith for (i=0;i<nz;i++) { 40db4efbfdSBarry Smith aj[i] = (unsigned short)AJ[i]; 41db4efbfdSBarry Smith } 42db4efbfdSBarry Smith inA->ops->setunfactored = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj; 43db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj; 44db4efbfdSBarry Smith ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr); 45db4efbfdSBarry Smith } else { 46db4efbfdSBarry Smith /* Scale the column indices for easier indexing in MatSolve. */ 47db4efbfdSBarry Smith /* for (i=0;i<nz;i++) { */ 48db4efbfdSBarry Smith /* AJ[i] = AJ[i]*4; */ 49db4efbfdSBarry Smith /* } */ 50db4efbfdSBarry Smith inA->ops->setunfactored = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE; 51db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE; 52db4efbfdSBarry Smith ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr); 53db4efbfdSBarry Smith } 54db4efbfdSBarry Smith # else 55db4efbfdSBarry Smith /* This should never be reached. If so, problem in PetscSSEIsEnabled. */ 56db4efbfdSBarry Smith SETERRQ(PETSC_ERR_SUP,"SSE Hardware unavailable"); 57db4efbfdSBarry Smith # endif 58db4efbfdSBarry Smith } else { 59db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 60db4efbfdSBarry Smith } 61db4efbfdSBarry Smith } 62db4efbfdSBarry Smith #else 63db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 64db4efbfdSBarry Smith #endif 65db4efbfdSBarry Smith break; 66db4efbfdSBarry Smith case 5: 67db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 68db4efbfdSBarry Smith break; 69db4efbfdSBarry Smith case 6: 70db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 71db4efbfdSBarry Smith break; 72db4efbfdSBarry Smith case 7: 73db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 74db4efbfdSBarry Smith break; 75db4efbfdSBarry Smith default: 76db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 77db4efbfdSBarry Smith break; 78db4efbfdSBarry Smith } 79db4efbfdSBarry Smith } else { 80db4efbfdSBarry Smith switch (inA->rmap->bs) { 81db4efbfdSBarry Smith case 1: 82db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 83db4efbfdSBarry Smith break; 84db4efbfdSBarry Smith case 2: 85db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 86db4efbfdSBarry Smith break; 87db4efbfdSBarry Smith case 3: 88db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 89db4efbfdSBarry Smith break; 90db4efbfdSBarry Smith case 4: 91db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 92db4efbfdSBarry Smith break; 93db4efbfdSBarry Smith case 5: 94db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 95db4efbfdSBarry Smith break; 96db4efbfdSBarry Smith case 6: 97db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 98db4efbfdSBarry Smith break; 99db4efbfdSBarry Smith case 7: 100db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7; 101db4efbfdSBarry Smith break; 102db4efbfdSBarry Smith default: 103db4efbfdSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 104db4efbfdSBarry Smith break; 105db4efbfdSBarry Smith } 106db4efbfdSBarry Smith } 107db4efbfdSBarry Smith PetscFunctionReturn(0); 108db4efbfdSBarry Smith } 109db4efbfdSBarry Smith 11083287d42SBarry Smith /* 11183287d42SBarry Smith The symbolic factorization code is identical to that for AIJ format, 11283287d42SBarry Smith except for very small changes since this is now a SeqBAIJ datastructure. 11383287d42SBarry Smith NOT good code reuse. 11483287d42SBarry Smith */ 115c7272abeSHong Zhang #include "petscbt.h" 116c7272abeSHong Zhang #include "../src/mat/utils/freespace.h" 117c7272abeSHong Zhang 1184a2ae208SSatish Balay #undef __FUNCT__ 119faca2338SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_newdatastruct" 120faca2338SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_newdatastruct(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 12183287d42SBarry Smith { 12283287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 123c7272abeSHong Zhang PetscInt n=a->mbs,bs = A->rmap->bs,bs2=a->bs2; 124c7272abeSHong Zhang PetscTruth row_identity,col_identity,both_identity; 12583287d42SBarry Smith IS isicol; 1266849ba73SBarry Smith PetscErrorCode ierr; 127c7272abeSHong Zhang const PetscInt *r,*ic; 128c7272abeSHong Zhang PetscInt i,*ai=a->i,*aj=a->j; 129c7272abeSHong Zhang PetscInt *bi,*bj,*ajtmp; 130c7272abeSHong Zhang PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 131c7272abeSHong Zhang PetscReal f; 132c7272abeSHong Zhang PetscInt nlnk,*lnk,k,**bi_ptr; 133c7272abeSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 134c7272abeSHong Zhang PetscBT lnkbt; 13583287d42SBarry Smith 13683287d42SBarry Smith PetscFunctionBegin; 137d0f46423SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 13883287d42SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 13983287d42SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 14083287d42SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 14183287d42SBarry Smith 14283287d42SBarry Smith /* get new row pointers */ 143faca2338SShri Abhyankar ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 144faca2338SShri Abhyankar bi[0] = 0; 145faca2338SShri Abhyankar 146faca2338SShri Abhyankar /* bdiag is location of diagonal in factor */ 147faca2338SShri Abhyankar ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 148faca2338SShri Abhyankar bdiag[0] = 0; 149faca2338SShri Abhyankar 150faca2338SShri Abhyankar /* linked list for storing column indices of the active row */ 151faca2338SShri Abhyankar nlnk = n + 1; 152faca2338SShri Abhyankar ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 153faca2338SShri Abhyankar 154faca2338SShri Abhyankar ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 155faca2338SShri Abhyankar 156faca2338SShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 157faca2338SShri Abhyankar f = info->fill; 158faca2338SShri Abhyankar ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 159faca2338SShri Abhyankar current_space = free_space; 160faca2338SShri Abhyankar 161faca2338SShri Abhyankar for (i=0; i<n; i++) { 162faca2338SShri Abhyankar /* copy previous fill into linked list */ 163faca2338SShri Abhyankar nzi = 0; 164faca2338SShri Abhyankar nnz = ai[r[i]+1] - ai[r[i]]; 165faca2338SShri Abhyankar if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 166faca2338SShri Abhyankar ajtmp = aj + ai[r[i]]; 167faca2338SShri Abhyankar ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 168faca2338SShri Abhyankar nzi += nlnk; 169faca2338SShri Abhyankar 170faca2338SShri Abhyankar /* add pivot rows into linked list */ 171faca2338SShri Abhyankar row = lnk[n]; 172faca2338SShri Abhyankar while (row < i) { 1738bfb2e2bSShri Abhyankar nzbd = bdiag[row] + 1; /* num of entries in the row with column index <= row */ 174faca2338SShri Abhyankar ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 175faca2338SShri Abhyankar ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 176faca2338SShri Abhyankar nzi += nlnk; 177faca2338SShri Abhyankar row = lnk[row]; 178faca2338SShri Abhyankar } 179faca2338SShri Abhyankar bi[i+1] = bi[i] + nzi; 180faca2338SShri Abhyankar im[i] = nzi; 181faca2338SShri Abhyankar 182faca2338SShri Abhyankar /* mark bdiag */ 183faca2338SShri Abhyankar nzbd = 0; 184faca2338SShri Abhyankar nnz = nzi; 185faca2338SShri Abhyankar k = lnk[n]; 186faca2338SShri Abhyankar while (nnz-- && k < i){ 187faca2338SShri Abhyankar nzbd++; 188faca2338SShri Abhyankar k = lnk[k]; 189faca2338SShri Abhyankar } 190783ef271SHong Zhang bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */ 191faca2338SShri Abhyankar 192faca2338SShri Abhyankar /* if free space is not available, make more free space */ 193faca2338SShri Abhyankar if (current_space->local_remaining<nzi) { 194faca2338SShri Abhyankar nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */ 195faca2338SShri Abhyankar ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 196faca2338SShri Abhyankar reallocs++; 197faca2338SShri Abhyankar } 198faca2338SShri Abhyankar 199faca2338SShri Abhyankar /* copy data into free space, then initialize lnk */ 200faca2338SShri Abhyankar ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 201faca2338SShri Abhyankar bi_ptr[i] = current_space->array; 202faca2338SShri Abhyankar current_space->array += nzi; 203faca2338SShri Abhyankar current_space->local_used += nzi; 204faca2338SShri Abhyankar current_space->local_remaining -= nzi; 205faca2338SShri Abhyankar } 206faca2338SShri Abhyankar #if defined(PETSC_USE_INFO) 207faca2338SShri Abhyankar if (ai[n] != 0) { 208faca2338SShri Abhyankar PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 209faca2338SShri Abhyankar ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 210faca2338SShri Abhyankar ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 211faca2338SShri Abhyankar ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 212faca2338SShri Abhyankar ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 213faca2338SShri Abhyankar } else { 214faca2338SShri Abhyankar ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 215faca2338SShri Abhyankar } 216faca2338SShri Abhyankar #endif 217faca2338SShri Abhyankar 218faca2338SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 219faca2338SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 220faca2338SShri Abhyankar 221faca2338SShri Abhyankar /* destroy list of free space and other temporary array(s) */ 222faca2338SShri Abhyankar ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 223783ef271SHong Zhang ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 224faca2338SShri Abhyankar ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 225faca2338SShri Abhyankar ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 226faca2338SShri Abhyankar 227faca2338SShri Abhyankar /* put together the new matrix */ 228faca2338SShri Abhyankar ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 229faca2338SShri Abhyankar ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 230faca2338SShri Abhyankar b = (Mat_SeqBAIJ*)(B)->data; 231faca2338SShri Abhyankar b->free_a = PETSC_TRUE; 232faca2338SShri Abhyankar b->free_ij = PETSC_TRUE; 233faca2338SShri Abhyankar b->singlemalloc = PETSC_FALSE; 234faca2338SShri Abhyankar ierr = PetscMalloc((bi[2*n+1])*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 235faca2338SShri Abhyankar b->j = bj; 236faca2338SShri Abhyankar b->i = bi; 237faca2338SShri Abhyankar b->diag = bdiag; 2387f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 239faca2338SShri Abhyankar b->ilen = 0; 240faca2338SShri Abhyankar b->imax = 0; 241faca2338SShri Abhyankar b->row = isrow; 242faca2338SShri Abhyankar b->col = iscol; 243faca2338SShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 244faca2338SShri Abhyankar ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 245faca2338SShri Abhyankar ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 246faca2338SShri Abhyankar b->icol = isicol; 247faca2338SShri Abhyankar ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 248faca2338SShri Abhyankar ierr = PetscLogObjectMemory(B,(bi[2*n+1])*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr); 249faca2338SShri Abhyankar 250faca2338SShri Abhyankar b->maxnz = b->nz = bi[2*n+1]; 2514c000e2eSHong Zhang B->factor = MAT_FACTOR_LU; 2524c000e2eSHong Zhang B->info.factor_mallocs = reallocs; 2534c000e2eSHong Zhang B->info.fill_ratio_given = f; 254faca2338SShri Abhyankar 255faca2338SShri Abhyankar if (ai[n] != 0) { 2564c000e2eSHong Zhang B->info.fill_ratio_needed = ((PetscReal)bi[2*n+1])/((PetscReal)ai[n]); 257faca2338SShri Abhyankar } else { 2584c000e2eSHong Zhang B->info.fill_ratio_needed = 0.0; 259faca2338SShri Abhyankar } 260b588c5a2SHong Zhang 261faca2338SShri Abhyankar ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 262faca2338SShri Abhyankar ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 263faca2338SShri Abhyankar both_identity = (PetscTruth) (row_identity && col_identity); 2648bfb2e2bSShri Abhyankar if(both_identity){ 2658bfb2e2bSShri Abhyankar switch(bs){ 2668bfb2e2bSShri Abhyankar case 2: 2674c000e2eSHong Zhang B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct; 2684c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering_newdatastruct; 2698bfb2e2bSShri Abhyankar break; 2708bfb2e2bSShri Abhyankar case 3: 271359bf53fSShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct; 2724c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_newdatastruct; 2738bfb2e2bSShri Abhyankar break; 2748bfb2e2bSShri Abhyankar case 4: 275209027a4SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_newdatastruct; 2764c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering_newdatastruct; 2778bfb2e2bSShri Abhyankar break; 2788bfb2e2bSShri Abhyankar case 5: 2790deeaf61SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_newdatastruct; 2804c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering_newdatastruct; 2818bfb2e2bSShri Abhyankar break; 2828bfb2e2bSShri Abhyankar case 6: 283*bef36659SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_newdatastruct; 2844c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering_newdatastruct; 2858bfb2e2bSShri Abhyankar break; 2868bfb2e2bSShri Abhyankar case 7: 287*bef36659SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_newdatastruct; 2884c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering_newdatastruct; 2898bfb2e2bSShri Abhyankar break; 2908bfb2e2bSShri Abhyankar default: 291*bef36659SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_newdatastruct; 2924c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct; 2938bfb2e2bSShri Abhyankar } 2948bfb2e2bSShri Abhyankar } 2958bfb2e2bSShri Abhyankar else { 2968bfb2e2bSShri Abhyankar switch(bs){ 2978bfb2e2bSShri Abhyankar case 2: 2984c000e2eSHong Zhang B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_newdatastruct; 2994c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_2_newdatastruct; 3008bfb2e2bSShri Abhyankar break; 3018bfb2e2bSShri Abhyankar case 3: 30217542729SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_newdatastruct; 3034c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_3_newdatastruct; 3048bfb2e2bSShri Abhyankar break; 3058bfb2e2bSShri Abhyankar case 4: 306209027a4SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_newdatastruct; 3074c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_4_newdatastruct; 3088bfb2e2bSShri Abhyankar break; 3098bfb2e2bSShri Abhyankar case 5: 3100deeaf61SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_newdatastruct; 3114c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_5_newdatastruct; 3128bfb2e2bSShri Abhyankar break; 3138bfb2e2bSShri Abhyankar case 6: 314*bef36659SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_newdatastruct; 3154c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_6_newdatastruct; 3168bfb2e2bSShri Abhyankar break; 3178bfb2e2bSShri Abhyankar case 7: 318*bef36659SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_newdatastruct; 3194c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_7_newdatastruct; 3208bfb2e2bSShri Abhyankar break; 3218bfb2e2bSShri Abhyankar default: 322*bef36659SShri Abhyankar B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_newdatastruct; 3234c000e2eSHong Zhang B->ops->solve = MatSolve_SeqBAIJ_N_newdatastruct; 3248bfb2e2bSShri Abhyankar } 3258bfb2e2bSShri Abhyankar } 326b588c5a2SHong Zhang /* ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr); */ 327faca2338SShri Abhyankar PetscFunctionReturn(0); 328faca2338SShri Abhyankar } 329faca2338SShri Abhyankar 330faca2338SShri Abhyankar #undef __FUNCT__ 331faca2338SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ" 332faca2338SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 333faca2338SShri Abhyankar { 334faca2338SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 335faca2338SShri Abhyankar PetscInt n=a->mbs,bs = A->rmap->bs,bs2=a->bs2; 336faca2338SShri Abhyankar PetscTruth row_identity,col_identity,both_identity; 337faca2338SShri Abhyankar IS isicol; 338faca2338SShri Abhyankar PetscErrorCode ierr; 339faca2338SShri Abhyankar const PetscInt *r,*ic; 340faca2338SShri Abhyankar PetscInt i,*ai=a->i,*aj=a->j; 341faca2338SShri Abhyankar PetscInt *bi,*bj,*ajtmp; 342faca2338SShri Abhyankar PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 343faca2338SShri Abhyankar PetscReal f; 344faca2338SShri Abhyankar PetscInt nlnk,*lnk,k,**bi_ptr; 345faca2338SShri Abhyankar PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 346faca2338SShri Abhyankar PetscBT lnkbt; 347faca2338SShri Abhyankar PetscTruth newdatastruct=PETSC_FALSE; 348faca2338SShri Abhyankar 349faca2338SShri Abhyankar PetscFunctionBegin; 350faca2338SShri Abhyankar ierr = PetscOptionsGetTruth(PETSC_NULL,"-lu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr); 351faca2338SShri Abhyankar if(newdatastruct){ 352faca2338SShri Abhyankar ierr = MatLUFactorSymbolic_SeqBAIJ_newdatastruct(B,A,isrow,iscol,info);CHKERRQ(ierr); 353faca2338SShri Abhyankar PetscFunctionReturn(0); 354faca2338SShri Abhyankar } 355faca2338SShri Abhyankar 356faca2338SShri Abhyankar if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 357faca2338SShri Abhyankar ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 358faca2338SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 359faca2338SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 360faca2338SShri Abhyankar 361faca2338SShri Abhyankar /* get new row pointers */ 362c7272abeSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 363c7272abeSHong Zhang bi[0] = 0; 364c7272abeSHong Zhang 365c7272abeSHong Zhang /* bdiag is location of diagonal in factor */ 366c7272abeSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 367c7272abeSHong Zhang bdiag[0] = 0; 368c7272abeSHong Zhang 369c7272abeSHong Zhang /* linked list for storing column indices of the active row */ 370c7272abeSHong Zhang nlnk = n + 1; 371c7272abeSHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 372c7272abeSHong Zhang 373c7272abeSHong Zhang ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 374c7272abeSHong Zhang 375c7272abeSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 376c7272abeSHong Zhang f = info->fill; 377c7272abeSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 378c7272abeSHong Zhang current_space = free_space; 37983287d42SBarry Smith 38083287d42SBarry Smith for (i=0; i<n; i++) { 381c7272abeSHong Zhang /* copy previous fill into linked list */ 38283287d42SBarry Smith nzi = 0; 383c7272abeSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 384c7272abeSHong Zhang if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 385c7272abeSHong Zhang ajtmp = aj + ai[r[i]]; 386c7272abeSHong Zhang ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 387c7272abeSHong Zhang nzi += nlnk; 388c7272abeSHong Zhang 389c7272abeSHong Zhang /* add pivot rows into linked list */ 390c7272abeSHong Zhang row = lnk[n]; 391c7272abeSHong Zhang while (row < i) { 392c7272abeSHong Zhang nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 393c7272abeSHong Zhang ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 394c7272abeSHong Zhang ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 395c7272abeSHong Zhang nzi += nlnk; 396c7272abeSHong Zhang row = lnk[row]; 39783287d42SBarry Smith } 398c7272abeSHong Zhang bi[i+1] = bi[i] + nzi; 399c7272abeSHong Zhang im[i] = nzi; 400c7272abeSHong Zhang 401c7272abeSHong Zhang /* mark bdiag */ 402c7272abeSHong Zhang nzbd = 0; 403c7272abeSHong Zhang nnz = nzi; 404c7272abeSHong Zhang k = lnk[n]; 405c7272abeSHong Zhang while (nnz-- && k < i){ 406c7272abeSHong Zhang nzbd++; 407c7272abeSHong Zhang k = lnk[k]; 408c7272abeSHong Zhang } 409c7272abeSHong Zhang bdiag[i] = bi[i] + nzbd; 410c7272abeSHong Zhang 411c7272abeSHong Zhang /* if free space is not available, make more free space */ 412c7272abeSHong Zhang if (current_space->local_remaining<nzi) { 413c7272abeSHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 414c7272abeSHong Zhang ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 415c7272abeSHong Zhang reallocs++; 41683287d42SBarry Smith } 41783287d42SBarry Smith 418c7272abeSHong Zhang /* copy data into free space, then initialize lnk */ 419c7272abeSHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 420c7272abeSHong Zhang bi_ptr[i] = current_space->array; 421c7272abeSHong Zhang current_space->array += nzi; 422c7272abeSHong Zhang current_space->local_used += nzi; 423c7272abeSHong Zhang current_space->local_remaining -= nzi; 424c7272abeSHong Zhang } 4256cf91177SBarry Smith #if defined(PETSC_USE_INFO) 42683287d42SBarry Smith if (ai[n] != 0) { 427c7272abeSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 428ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 429ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 430ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 431ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 43283287d42SBarry Smith } else { 433c7272abeSHong Zhang ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 43483287d42SBarry Smith } 43563ba0a88SBarry Smith #endif 43683287d42SBarry Smith 43783287d42SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 43883287d42SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 43983287d42SBarry Smith 440c7272abeSHong Zhang /* destroy list of free space and other temporary array(s) */ 441c7272abeSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 442c7272abeSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 443c7272abeSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 444c7272abeSHong Zhang ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 44583287d42SBarry Smith 44683287d42SBarry Smith /* put together the new matrix */ 447719d5645SBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 448719d5645SBarry Smith ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 449719d5645SBarry Smith b = (Mat_SeqBAIJ*)(B)->data; 450e6b907acSBarry Smith b->free_a = PETSC_TRUE; 451e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 452c7272abeSHong Zhang b->singlemalloc = PETSC_FALSE; 453c7272abeSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 454c7272abeSHong Zhang b->j = bj; 455c7272abeSHong Zhang b->i = bi; 456c7272abeSHong Zhang b->diag = bdiag; 4577f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 45883287d42SBarry Smith b->ilen = 0; 45983287d42SBarry Smith b->imax = 0; 46083287d42SBarry Smith b->row = isrow; 46183287d42SBarry Smith b->col = iscol; 462d3d32019SBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 46383287d42SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 46483287d42SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 46583287d42SBarry Smith b->icol = isicol; 46687828ca2SBarry Smith ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 467c7272abeSHong Zhang ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr); 46883287d42SBarry Smith 469c7272abeSHong Zhang b->maxnz = b->nz = bi[n] ; 470c7272abeSHong Zhang (B)->factor = MAT_FACTOR_LU; 471719d5645SBarry Smith (B)->info.factor_mallocs = reallocs; 472719d5645SBarry Smith (B)->info.fill_ratio_given = f; 473c7272abeSHong Zhang 47483287d42SBarry Smith if (ai[n] != 0) { 475c7272abeSHong Zhang (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 47683287d42SBarry Smith } else { 477719d5645SBarry Smith (B)->info.fill_ratio_needed = 0.0; 47883287d42SBarry Smith } 479c7272abeSHong Zhang 480db4efbfdSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 481db4efbfdSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 4827d18ce8fSMatthew Knepley both_identity = (PetscTruth) (row_identity && col_identity); 48341df41f0SMatthew Knepley ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr); 48483287d42SBarry Smith PetscFunctionReturn(0); 48583287d42SBarry Smith } 486db4efbfdSBarry Smith 487