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