xref: /petsc/src/mat/impls/aij/seq/bas/basfactor.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1b2f1ab58SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h>
5b2f1ab58SBarry Smith 
6b2f1ab58SBarry Smith #undef __FUNCT__
7b2f1ab58SBarry Smith #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ_Bas"
8b2f1ab58SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
9b2f1ab58SBarry Smith {
10b2f1ab58SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
11b2f1ab58SBarry Smith   Mat_SeqSBAIJ   *b;
12b2f1ab58SBarry Smith   PetscErrorCode ierr;
13ace3abfcSBarry Smith   PetscBool      perm_identity,missing;
14b2f1ab58SBarry Smith   PetscInt       reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
15b2f1ab58SBarry Smith   const PetscInt *rip,*riip;
16b2f1ab58SBarry Smith   PetscInt       j;
17b2f1ab58SBarry Smith   PetscInt       d;
18b2f1ab58SBarry Smith   PetscInt       ncols,*cols,*uj;
19b2f1ab58SBarry Smith   PetscReal      fill=info->fill,levels=info->levels;
20b2f1ab58SBarry Smith   IS             iperm;
21b2f1ab58SBarry Smith   spbas_matrix   Pattern_0, Pattern_P;
22b2f1ab58SBarry Smith 
23b2f1ab58SBarry Smith   PetscFunctionBegin;
24e32f2f54SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
25b2f1ab58SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
26e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
27b2f1ab58SBarry Smith   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
28b2f1ab58SBarry Smith   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
29b2f1ab58SBarry Smith 
30b2f1ab58SBarry Smith   /* ICC(0) without matrix ordering: simply copies fill pattern */
31b2f1ab58SBarry Smith   if (!levels && perm_identity) {
32b2f1ab58SBarry Smith     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
33b2f1ab58SBarry Smith     ui[0] = 0;
34b2f1ab58SBarry Smith 
35b2f1ab58SBarry Smith     for (i=0; i<am; i++) {
36b2f1ab58SBarry Smith       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
37b2f1ab58SBarry Smith     }
38b2f1ab58SBarry Smith     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
39b2f1ab58SBarry Smith     cols = uj;
40b2f1ab58SBarry Smith     for (i=0; i<am; i++) {
41b2f1ab58SBarry Smith       aj    = a->j + a->diag[i];
42b2f1ab58SBarry Smith       ncols = ui[i+1] - ui[i];
43b2f1ab58SBarry Smith       for (j=0; j<ncols; j++) *cols++ = *aj++;
44b2f1ab58SBarry Smith     }
45b2f1ab58SBarry Smith   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
46b2f1ab58SBarry Smith     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
47b2f1ab58SBarry Smith     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
48b2f1ab58SBarry Smith 
4971d55d18SBarry Smith     /* Create spbas_matrix for pattern */
50b2f1ab58SBarry Smith     ierr = spbas_pattern_only(am, am, ai, aj, &Pattern_0);CHKERRQ(ierr);
51b2f1ab58SBarry Smith 
5271d55d18SBarry Smith     /* Apply the permutation */
53b2f1ab58SBarry Smith     ierr = spbas_apply_reordering(&Pattern_0, rip, riip);CHKERRQ(ierr);
54b2f1ab58SBarry Smith 
5571d55d18SBarry Smith     /* Raise the power */
5671d55d18SBarry Smith     ierr = spbas_power(Pattern_0, (int) levels+1, &Pattern_P);CHKERRQ(ierr);
57b2f1ab58SBarry Smith     ierr = spbas_delete(Pattern_0);CHKERRQ(ierr);
58b2f1ab58SBarry Smith 
5971d55d18SBarry Smith     /* Keep only upper triangle of pattern */
6022d28d08SBarry Smith     ierr = spbas_keep_upper(&Pattern_P);CHKERRQ(ierr);
61b2f1ab58SBarry Smith 
6271d55d18SBarry Smith     /* Convert to Sparse Row Storage  */
63c328eeadSBarry Smith     ierr = spbas_matrix_to_crs(Pattern_P, PETSC_NULL, &ui, &uj);CHKERRQ(ierr);
64b2f1ab58SBarry Smith     ierr = spbas_delete(Pattern_P);CHKERRQ(ierr);
65b2f1ab58SBarry Smith   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
66b2f1ab58SBarry Smith 
67b2f1ab58SBarry Smith   /* put together the new matrix in MATSEQSBAIJ format */
68b2f1ab58SBarry Smith 
69b2f1ab58SBarry Smith   b                = (Mat_SeqSBAIJ*)(fact)->data;
70b2f1ab58SBarry Smith   b->singlemalloc  = PETSC_FALSE;
71*2205254eSKarl Rupp 
72b2f1ab58SBarry Smith   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
73*2205254eSKarl Rupp 
74b2f1ab58SBarry Smith   b->j             = uj;
75b2f1ab58SBarry Smith   b->i             = ui;
76b2f1ab58SBarry Smith   b->diag          = 0;
77b2f1ab58SBarry Smith   b->ilen          = 0;
78b2f1ab58SBarry Smith   b->imax          = 0;
79b2f1ab58SBarry Smith   b->row           = perm;
80b2f1ab58SBarry Smith   b->col           = perm;
81*2205254eSKarl Rupp 
82b2f1ab58SBarry Smith   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
83b2f1ab58SBarry Smith   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
84*2205254eSKarl Rupp 
85b2f1ab58SBarry Smith   b->icol          = iperm;
86b2f1ab58SBarry Smith   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
87b2f1ab58SBarry Smith   ierr             = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
88b2f1ab58SBarry Smith   ierr             = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
89b2f1ab58SBarry Smith   b->maxnz         = b->nz = ui[am];
90b2f1ab58SBarry Smith   b->free_a        = PETSC_TRUE;
91b2f1ab58SBarry Smith   b->free_ij       = PETSC_TRUE;
92b2f1ab58SBarry Smith 
93b2f1ab58SBarry Smith   (fact)->info.factor_mallocs   = reallocs;
94b2f1ab58SBarry Smith   (fact)->info.fill_ratio_given = fill;
95b2f1ab58SBarry Smith   if (ai[am] != 0) {
96b2f1ab58SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
97b2f1ab58SBarry Smith   } else {
98b2f1ab58SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
99b2f1ab58SBarry Smith   }
100b2f1ab58SBarry Smith   /*  (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */
101b2f1ab58SBarry Smith   PetscFunctionReturn(0);
102b2f1ab58SBarry Smith }
103b2f1ab58SBarry Smith 
104b2f1ab58SBarry Smith 
105b2f1ab58SBarry Smith #undef __FUNCT__
106b2f1ab58SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_Bas"
107b2f1ab58SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo *info)
108b2f1ab58SBarry Smith {
109b2f1ab58SBarry Smith   Mat            C = B;
110b2f1ab58SBarry Smith   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
111b2f1ab58SBarry Smith   IS             ip=b->row,iip = b->icol;
112b2f1ab58SBarry Smith   PetscErrorCode ierr;
113b2f1ab58SBarry Smith   const PetscInt *rip,*riip;
114b2f1ab58SBarry Smith   PetscInt       mbs=A->rmap->n,*bi=b->i,*bj=b->j;
115b2f1ab58SBarry Smith 
116b2f1ab58SBarry Smith   MatScalar    *ba     = b->a;
117c88783f4SHong Zhang   PetscReal    shiftnz = info->shiftamount;
118d36aa34eSBarry Smith   PetscReal    droptol = -1;
119ace3abfcSBarry Smith   PetscBool    perm_identity;
120b2f1ab58SBarry Smith   spbas_matrix Pattern, matrix_L,matrix_LT;
1219767453cSBarry Smith   PetscReal    mem_reduction;
122b2f1ab58SBarry Smith 
123b2f1ab58SBarry Smith   PetscFunctionBegin;
12471d55d18SBarry Smith   /* Reduce memory requirements:   erase values of B-matrix */
125b2f1ab58SBarry Smith   ierr = PetscFree(ba);CHKERRQ(ierr);
12671d55d18SBarry Smith   /*   Compress (maximum) sparseness pattern of B-matrix */
1279767453cSBarry Smith   ierr = spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);CHKERRQ(ierr);
128b2f1ab58SBarry Smith   ierr = PetscFree(bi);CHKERRQ(ierr);
129b2f1ab58SBarry Smith   ierr = PetscFree(bj);CHKERRQ(ierr);
130b2f1ab58SBarry Smith 
1319767453cSBarry Smith   ierr = PetscInfo1(PETSC_NULL,"    compression rate for spbas_compress_pattern %G \n",mem_reduction);CHKERRQ(ierr);
132b2f1ab58SBarry Smith 
1334e1ff37aSBarry Smith   /* Make Cholesky decompositions with larger Manteuffel shifts until no more    negative diagonals are found. */
134b2f1ab58SBarry Smith   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
135b2f1ab58SBarry Smith   ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr);
136b2f1ab58SBarry Smith 
137b2f1ab58SBarry Smith   if (info->usedt) {
138b2f1ab58SBarry Smith     droptol = info->dt;
139b2f1ab58SBarry Smith   }
140b2f1ab58SBarry Smith   for (ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;)
141b2f1ab58SBarry Smith   {
14271d55d18SBarry Smith     ierr = spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);CHKERRQ(ierr);
143*2205254eSKarl Rupp     if (ierr == NEGATIVE_DIAGONAL) {
144b2f1ab58SBarry Smith       shiftnz *= 1.5;
1459767453cSBarry Smith       if (shiftnz < 1e-5) shiftnz=1e-5;
1469767453cSBarry Smith       ierr = PetscInfo1(PETSC_NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%G\n",shiftnz);CHKERRQ(ierr);
147b2f1ab58SBarry Smith     }
148b2f1ab58SBarry Smith   }
149b2f1ab58SBarry Smith   ierr = spbas_delete(Pattern);CHKERRQ(ierr);
150b2f1ab58SBarry Smith 
1514e1ff37aSBarry Smith   ierr = PetscInfo1(PETSC_NULL,"    memory_usage for  spbas_incomplete_cholesky  %G bytes per row\n", (PetscReal) spbas_memory_requirement(matrix_LT)/ (PetscReal) mbs);CHKERRQ(ierr);
152b2f1ab58SBarry Smith 
153b2f1ab58SBarry Smith   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
154b2f1ab58SBarry Smith   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
155b2f1ab58SBarry Smith 
15671d55d18SBarry Smith   /* Convert spbas_matrix to compressed row storage */
157b2f1ab58SBarry Smith   ierr = spbas_transpose(matrix_LT, &matrix_L);CHKERRQ(ierr);
158b2f1ab58SBarry Smith   ierr = spbas_delete(matrix_LT);CHKERRQ(ierr);
159b2f1ab58SBarry Smith   ierr = spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);CHKERRQ(ierr);
160b2f1ab58SBarry Smith   b->i =bi; b->j=bj; b->a=ba;
161b2f1ab58SBarry Smith   ierr = spbas_delete(matrix_L);CHKERRQ(ierr);
162b2f1ab58SBarry Smith 
16371d55d18SBarry Smith   /* Set the appropriate solution functions */
164b2f1ab58SBarry Smith   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
165b2f1ab58SBarry Smith   if (perm_identity) {
166b2f1ab58SBarry Smith     (B)->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
167b2f1ab58SBarry Smith     (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
168b2f1ab58SBarry Smith     (B)->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
169b2f1ab58SBarry Smith     (B)->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
170b2f1ab58SBarry Smith   } else {
171b2f1ab58SBarry Smith     (B)->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
172b2f1ab58SBarry Smith     (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
173b2f1ab58SBarry Smith     (B)->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
174b2f1ab58SBarry Smith     (B)->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
175b2f1ab58SBarry Smith   }
176b2f1ab58SBarry Smith 
177b2f1ab58SBarry Smith   C->assembled    = PETSC_TRUE;
178b2f1ab58SBarry Smith   C->preallocated = PETSC_TRUE;
179*2205254eSKarl Rupp 
180b2f1ab58SBarry Smith   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
181b2f1ab58SBarry Smith   PetscFunctionReturn(0);
182b2f1ab58SBarry Smith }
183b2f1ab58SBarry Smith 
1841e364ec1SMatthew G Knepley EXTERN_C_BEGIN
185b2f1ab58SBarry Smith #undef __FUNCT__
186b2f1ab58SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_bas"
187b2f1ab58SBarry Smith PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B)
188b2f1ab58SBarry Smith {
189b2f1ab58SBarry Smith   PetscInt       n = A->rmap->n;
190b2f1ab58SBarry Smith   PetscErrorCode ierr;
191b2f1ab58SBarry Smith 
192b2f1ab58SBarry Smith   PetscFunctionBegin;
193b2f1ab58SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
194b2f1ab58SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
195b2f1ab58SBarry Smith   if (ftype == MAT_FACTOR_ICC) {
196b2f1ab58SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
197b2f1ab58SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
198*2205254eSKarl Rupp 
199b2f1ab58SBarry Smith     (*B)->ops->iccfactorsymbolic     = MatICCFactorSymbolic_SeqAIJ_Bas;
200b2f1ab58SBarry Smith     (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
201e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
202d5f3da31SBarry Smith   (*B)->factortype = ftype;
203b2f1ab58SBarry Smith   PetscFunctionReturn(0);
204b2f1ab58SBarry Smith }
2051e364ec1SMatthew G Knepley EXTERN_C_END
206