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