xref: /petsc/src/mat/impls/aij/seq/bas/basfactor.c (revision 71d55d18a6cc105fecd21255dfee54cec910a2af)
1 
2 
3 #include <petsc.h>
4 #include <petscksp.h>
5 #include "private/kspimpl.h"
6 #include "petscpc.h"
7 #include "../src/mat/impls/aij/seq/aij.h"
8 #include "../src/mat/impls/sbaij/seq/sbaij.h"
9 #include "../src/mat/impls/aij/seq/bas/spbas.h"
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ_Bas"
13 PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
14 {
15   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
16   Mat_SeqSBAIJ       *b;
17   PetscErrorCode     ierr;
18   PetscTruth         perm_identity,missing;
19   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
20   const PetscInt     *rip,*riip;
21   PetscInt           j;
22   PetscInt           d;
23   PetscInt           ncols,*cols,*uj;
24   PetscReal          fill=info->fill,levels=info->levels;
25   IS                 iperm;
26   spbas_matrix       Pattern_0, Pattern_P;
27 
28   PetscFunctionBegin;
29   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);
30   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
31   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
32   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
33   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
34 
35 
36   /* ICC(0) without matrix ordering: simply copies fill pattern */
37   if (!levels && perm_identity) {
38     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
39     ui[0] = 0;
40 
41     for (i=0; i<am; i++) {
42       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
43     }
44     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
45     cols = uj;
46     for (i=0; i<am; i++) {
47       aj    = a->j + a->diag[i];
48       ncols = ui[i+1] - ui[i];
49       for (j=0; j<ncols; j++) *cols++ = *aj++;
50     }
51   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
52     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
53     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
54 
55     /* Create spbas_matrix for pattern */
56     ierr = spbas_pattern_only(am, am, ai, aj, &Pattern_0);CHKERRQ(ierr);
57 
58     /* Apply the permutation */
59     ierr = spbas_apply_reordering( &Pattern_0, rip, riip);CHKERRQ(ierr);
60 
61     /* Raise the power */
62     ierr = spbas_power( Pattern_0, (int) levels+1, &Pattern_P);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, PETSC_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   PetscReal    mem_reduction;
124 
125   PetscFunctionBegin;
126   /* Reduce memory requirements:   erase values of B-matrix */
127   ierr = PetscFree(ba);CHKERRQ(ierr);
128   /*   Compress (maximum) sparseness pattern of B-matrix */
129   ierr = spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);CHKERRQ(ierr);
130   ierr = PetscFree(bi);CHKERRQ(ierr);
131   ierr = PetscFree(bj);CHKERRQ(ierr);
132 
133   ierr = PetscInfo1(PETSC_NULL,"    compression rate for spbas_compress_pattern %G \n",mem_reduction);CHKERRQ(ierr);
134 
135   /* Make Cholesky decompositions with larger Manteuffel shifts until no more
136      negative diagonals are found. */
137   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
138   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
139 
140   if (info->usedt) {
141     droptol = info->dt;
142   }
143   for (ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL; )
144   {
145      ierr  = spbas_incomplete_cholesky( A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);CHKERRQ(ierr);
146      if (ierr == NEGATIVE_DIAGONAL)
147      {
148         shiftnz *= 1.5;
149         if (shiftnz < 1e-5) shiftnz=1e-5;
150         ierr = PetscInfo1(PETSC_NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%G\n",shiftnz);CHKERRQ(ierr);
151      }
152   }
153   ierr = spbas_delete(Pattern);CHKERRQ(ierr);
154 
155   ierr = PetscInfo1(PETSC_NULL,"    memory_usage for  spbas_incomplete_cholesky  %G bytes per row\n",
156               (PetscReal) spbas_memory_requirement( matrix_LT)/ (PetscReal) mbs);CHKERRQ(ierr);
157 
158   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
159   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
160 
161   /* Convert spbas_matrix to compressed row storage */
162   ierr = spbas_transpose(matrix_LT, &matrix_L);CHKERRQ(ierr);
163   ierr = spbas_delete(matrix_LT);CHKERRQ(ierr);
164   ierr = spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);CHKERRQ(ierr);
165   b->i=bi; b->j=bj; b->a=ba;
166   ierr = spbas_delete(matrix_L);CHKERRQ(ierr);
167 
168   /* Set the appropriate solution functions */
169   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
170   if (perm_identity){
171     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
172     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
173     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
174     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
175   } else {
176     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_inplace;
177     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_inplace;
178     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_inplace;
179     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_inplace;
180   }
181 
182   C->assembled    = PETSC_TRUE;
183   C->preallocated = PETSC_TRUE;
184   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
185   PetscFunctionReturn(0);
186 }
187 
188 EXTERN_C_BEGIN
189 #undef __FUNCT__
190 #define __FUNCT__ "MatGetFactor_seqaij_bas"
191 PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B)
192 {
193   PetscInt           n = A->rmap->n;
194   PetscErrorCode     ierr;
195 
196   PetscFunctionBegin;
197   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
198   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
199   if (ftype == MAT_FACTOR_ICC) {
200     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
201     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
202     (*B)->ops->iccfactorsymbolic     = MatICCFactorSymbolic_SeqAIJ_Bas;
203     (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
204   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
205   (*B)->factor = ftype;
206   PetscFunctionReturn(0);
207 }
208 EXTERN_C_END
209