xref: /petsc/src/mat/impls/aij/seq/bas/basfactor.c (revision 3923b477fd0dced8a2d147b4fb4519fe3af97d3f)
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 = PetscMalloc((am+1)*sizeof(PetscInt),&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 = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&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 );
61 
62     /* Convert to Sparse Row Storage  */
63     ierr = spbas_matrix_to_crs(Pattern_P, PETSC_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   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
72   b->j    = uj;
73   b->i    = ui;
74   b->diag = 0;
75   b->ilen = 0;
76   b->imax = 0;
77   b->row  = perm;
78   b->col  = perm;
79   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
80   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
81   b->icol = iperm;
82   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
83   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
84   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
85   b->maxnz   = b->nz = ui[am];
86   b->free_a  = PETSC_TRUE;
87   b->free_ij = PETSC_TRUE;
88 
89   (fact)->info.factor_mallocs    = reallocs;
90   (fact)->info.fill_ratio_given  = fill;
91   if (ai[am] != 0) {
92     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
93   } else {
94     (fact)->info.fill_ratio_needed = 0.0;
95   }
96   /*  (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */
97   PetscFunctionReturn(0);
98 }
99 
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_Bas"
103 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo *info)
104 {
105   Mat            C = B;
106   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
107   IS             ip=b->row,iip = b->icol;
108   PetscErrorCode ierr;
109   const PetscInt *rip,*riip;
110   PetscInt       mbs=A->rmap->n,*bi=b->i,*bj=b->j;
111 
112   MatScalar      *ba=b->a;
113   PetscReal      shiftnz = info->shiftamount;
114   PetscReal      droptol = -1;
115   PetscBool      perm_identity;
116   spbas_matrix   Pattern, matrix_L,matrix_LT;
117   PetscReal      mem_reduction;
118 
119   PetscFunctionBegin;
120   /* Reduce memory requirements:   erase values of B-matrix */
121   ierr = PetscFree(ba);CHKERRQ(ierr);
122   /*   Compress (maximum) sparseness pattern of B-matrix */
123   ierr = spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);CHKERRQ(ierr);
124   ierr = PetscFree(bi);CHKERRQ(ierr);
125   ierr = PetscFree(bj);CHKERRQ(ierr);
126 
127   ierr = PetscInfo1(PETSC_NULL,"    compression rate for spbas_compress_pattern %G \n",mem_reduction);CHKERRQ(ierr);
128 
129   /* Make Cholesky decompositions with larger Manteuffel shifts until no more    negative diagonals are found. */
130   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
131   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
132 
133   if (info->usedt) {
134     droptol = info->dt;
135   }
136   for (ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL; )
137   {
138      ierr  = spbas_incomplete_cholesky( A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);CHKERRQ(ierr);
139      if (ierr == NEGATIVE_DIAGONAL)
140      {
141         shiftnz *= 1.5;
142         if (shiftnz < 1e-5) shiftnz=1e-5;
143         ierr = PetscInfo1(PETSC_NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%G\n",shiftnz);CHKERRQ(ierr);
144      }
145   }
146   ierr = spbas_delete(Pattern);CHKERRQ(ierr);
147 
148   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);
149 
150   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
151   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
152 
153   /* Convert spbas_matrix to compressed row storage */
154   ierr = spbas_transpose(matrix_LT, &matrix_L);CHKERRQ(ierr);
155   ierr = spbas_delete(matrix_LT);CHKERRQ(ierr);
156   ierr = spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);CHKERRQ(ierr);
157   b->i=bi; b->j=bj; b->a=ba;
158   ierr = spbas_delete(matrix_L);CHKERRQ(ierr);
159 
160   /* Set the appropriate solution functions */
161   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
162   if (perm_identity){
163     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
164     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
165     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
166     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
167   } else {
168     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_inplace;
169     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_inplace;
170     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_inplace;
171     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_inplace;
172   }
173 
174   C->assembled    = PETSC_TRUE;
175   C->preallocated = PETSC_TRUE;
176   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
177   PetscFunctionReturn(0);
178 }
179 
180 EXTERN_C_BEGIN
181 #undef __FUNCT__
182 #define __FUNCT__ "MatGetFactor_seqaij_bas"
183 PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B)
184 {
185   PetscInt           n = A->rmap->n;
186   PetscErrorCode     ierr;
187 
188   PetscFunctionBegin;
189   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
190   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
191   if (ftype == MAT_FACTOR_ICC) {
192     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
193     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
194     (*B)->ops->iccfactorsymbolic     = MatICCFactorSymbolic_SeqAIJ_Bas;
195     (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
196   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
197   (*B)->factortype = ftype;
198   PetscFunctionReturn(0);
199 }
200 EXTERN_C_END
201