Actual source code: basfactor.c
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
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: PetscTruth 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;
24: 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);
25: MatMissingDiagonal(A,&missing,&d);
26: if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
27: ISIdentity(perm,&perm_identity);
28: ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
30: /* ICC(0) without matrix ordering: simply copies fill pattern */
31: if (!levels && perm_identity) {
32: PetscMalloc((am+1)*sizeof(PetscInt),&ui);
33: ui[0] = 0;
35: for (i=0; i<am; i++) {
36: ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
37: }
38: PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);
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: ISGetIndices(iperm,&riip);
47: ISGetIndices(perm,&rip);
49: /* Create spbas_matrix for pattern */
50: spbas_pattern_only(am, am, ai, aj, &Pattern_0);
52: /* Apply the permutation */
53: spbas_apply_reordering( &Pattern_0, rip, riip);
54:
55: /* Raise the power */
56: spbas_power( Pattern_0, (int) levels+1, &Pattern_P);
57: spbas_delete( Pattern_0 );
59: /* Keep only upper triangle of pattern */
60: spbas_keep_upper( &Pattern_P );
62: /* Convert to Sparse Row Storage */
63: spbas_matrix_to_crs(Pattern_P, PETSC_NULL, &ui, &uj);
64: spbas_delete(Pattern_P);
65: } /* end of case: levels>0 || (levels=0 && !perm_identity) */
67: /* put together the new matrix in MATSEQSBAIJ format */
69: b = (Mat_SeqSBAIJ*)(fact)->data;
70: b->singlemalloc = PETSC_FALSE;
71: PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);
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: PetscObjectReference((PetscObject)perm);
80: PetscObjectReference((PetscObject)perm);
81: b->icol = iperm;
82: b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
83: PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);
84: PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
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: return(0);
98: }
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;
109: const PetscInt *rip,*riip;
110: PetscInt mbs=A->rmap->n,*bi=b->i,*bj=b->j;
112: MatScalar *ba=b->a;
113: PetscReal shiftnz = info->shiftamount;
114: PetscReal droptol = -1;
115: PetscTruth perm_identity;
116: spbas_matrix Pattern, matrix_L,matrix_LT;
117: PetscReal mem_reduction;
120: /* Reduce memory requirements: erase values of B-matrix */
121: PetscFree(ba);
122: /* Compress (maximum) sparseness pattern of B-matrix */
123: spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);
124: PetscFree(bi);
125: PetscFree(bj);
127: PetscInfo1(PETSC_NULL," compression rate for spbas_compress_pattern %G \n",mem_reduction);
129: /* Make Cholesky decompositions with larger Manteuffel shifts until no more negative diagonals are found. */
130: ISGetIndices(ip,&rip);
131: ISGetIndices(iip,&riip);
133: if (info->usedt) {
134: droptol = info->dt;
135: }
136: for (NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL; )
137: {
138: spbas_incomplete_cholesky( A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);
139: if (ierr == NEGATIVE_DIAGONAL)
140: {
141: shiftnz *= 1.5;
142: if (shiftnz < 1e-5) shiftnz=1e-5;
143: PetscInfo1(PETSC_NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%G\n",shiftnz);
144: }
145: }
146: spbas_delete(Pattern);
148: PetscInfo1(PETSC_NULL," memory_usage for spbas_incomplete_cholesky %G bytes per row\n", (PetscReal) spbas_memory_requirement( matrix_LT)/ (PetscReal) mbs);
150: ISRestoreIndices(ip,&rip);
151: ISRestoreIndices(iip,&riip);
153: /* Convert spbas_matrix to compressed row storage */
154: spbas_transpose(matrix_LT, &matrix_L);
155: spbas_delete(matrix_LT);
156: spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);
157: b->i=bi; b->j=bj; b->a=ba;
158: spbas_delete(matrix_L);
160: /* Set the appropriate solution functions */
161: ISIdentity(ip,&perm_identity);
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: }
174: C->assembled = PETSC_TRUE;
175: C->preallocated = PETSC_TRUE;
176: PetscLogFlops(C->rmap->n);
177: return(0);
178: }
183: PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B)
184: {
185: PetscInt n = A->rmap->n;
186: PetscErrorCode ierr;
189: MatCreate(((PetscObject)A)->comm,B);
190: MatSetSizes(*B,n,n,n,n);
191: if (ftype == MAT_FACTOR_ICC) {
192: MatSetType(*B,MATSEQSBAIJ);
193: MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);
194: (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas;
195: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
196: } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
197: (*B)->factor = ftype;
198: return(0);
199: }