Actual source code: basfactor.c
petsc-3.13.6 2020-09-29
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>
6: PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
7: {
8: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
9: Mat_SeqSBAIJ *b;
11: PetscBool perm_identity,missing;
12: PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
13: const PetscInt *rip,*riip;
14: PetscInt j;
15: PetscInt d;
16: PetscInt ncols,*cols,*uj;
17: PetscReal fill=info->fill,levels=info->levels;
18: IS iperm;
19: spbas_matrix Pattern_0, Pattern_P;
22: 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);
23: MatMissingDiagonal(A,&missing,&d);
24: if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
25: ISIdentity(perm,&perm_identity);
26: ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
28: /* ICC(0) without matrix ordering: simply copies fill pattern */
29: if (!levels && perm_identity) {
30: PetscMalloc1(am+1,&ui);
31: ui[0] = 0;
33: for (i=0; i<am; i++) {
34: ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
35: }
36: PetscMalloc1(ui[am]+1,&uj);
37: cols = uj;
38: for (i=0; i<am; i++) {
39: aj = a->j + a->diag[i];
40: ncols = ui[i+1] - ui[i];
41: for (j=0; j<ncols; j++) *cols++ = *aj++;
42: }
43: } else { /* case: levels>0 || (levels=0 && !perm_identity) */
44: ISGetIndices(iperm,&riip);
45: ISGetIndices(perm,&rip);
47: /* Create spbas_matrix for pattern */
48: spbas_pattern_only(am, am, ai, aj, &Pattern_0);
50: /* Apply the permutation */
51: spbas_apply_reordering(&Pattern_0, rip, riip);
53: /* Raise the power */
54: spbas_power(Pattern_0, (int) levels+1, &Pattern_P);
55: spbas_delete(Pattern_0);
57: /* Keep only upper triangle of pattern */
58: spbas_keep_upper(&Pattern_P);
60: /* Convert to Sparse Row Storage */
61: spbas_matrix_to_crs(Pattern_P, NULL, &ui, &uj);
62: spbas_delete(Pattern_P);
63: } /* end of case: levels>0 || (levels=0 && !perm_identity) */
65: /* put together the new matrix in MATSEQSBAIJ format */
67: b = (Mat_SeqSBAIJ*)(fact)->data;
68: b->singlemalloc = PETSC_FALSE;
70: PetscMalloc1(ui[am]+1,&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;
80: PetscObjectReference((PetscObject)perm);
81: PetscObjectReference((PetscObject)perm);
83: b->icol = iperm;
84: b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
85: PetscMalloc1(am+1,&b->solve_work);
86: PetscLogObjectMemory((PetscObject)(fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
87: b->maxnz = b->nz = ui[am];
88: b->free_a = PETSC_TRUE;
89: b->free_ij = PETSC_TRUE;
91: (fact)->info.factor_mallocs = reallocs;
92: (fact)->info.fill_ratio_given = fill;
93: if (ai[am] != 0) {
94: (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
95: } else {
96: (fact)->info.fill_ratio_needed = 0.0;
97: }
98: /* (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */
99: return(0);
100: }
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: PetscBool 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(NULL," compression rate for spbas_compress_pattern %g \n",(double)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: shiftnz *= 1.5;
141: if (shiftnz < 1e-5) shiftnz=1e-5;
142: PetscInfo1(NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n",(double)shiftnz);
143: }
144: }
145: spbas_delete(Pattern);
147: PetscInfo1(NULL," memory_usage for spbas_incomplete_cholesky %g bytes per row\n", (double)(PetscReal) (spbas_memory_requirement(matrix_LT)/ (PetscReal) mbs));
149: ISRestoreIndices(ip,&rip);
150: ISRestoreIndices(iip,&riip);
152: /* Convert spbas_matrix to compressed row storage */
153: spbas_transpose(matrix_LT, &matrix_L);
154: spbas_delete(matrix_LT);
155: spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);
156: b->i =bi; b->j=bj; b->a=ba;
157: spbas_delete(matrix_L);
159: /* Set the appropriate solution functions */
160: ISIdentity(ip,&perm_identity);
161: if (perm_identity) {
162: (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
163: (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
164: (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
165: (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
166: } else {
167: (B)->ops->solve = MatSolve_SeqSBAIJ_1_inplace;
168: (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
169: (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace;
170: (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace;
171: }
173: C->assembled = PETSC_TRUE;
174: C->preallocated = PETSC_TRUE;
176: PetscLogFlops(C->rmap->n);
177: return(0);
178: }
180: PetscErrorCode MatFactorGetSolverType_seqaij_bas(Mat A,MatSolverType *type)
181: {
183: *type = MATSOLVERBAS;
184: return(0);
185: }
187: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B)
188: {
189: PetscInt n = A->rmap->n;
193: MatCreate(PetscObjectComm((PetscObject)A),B);
194: MatSetSizes(*B,n,n,n,n);
195: if (ftype == MAT_FACTOR_ICC) {
196: MatSetType(*B,MATSEQSBAIJ);
197: MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,NULL);
199: (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas;
200: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
201: PetscObjectComposeFunction((PetscObject)*B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_bas);
202: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
203: (*B)->factortype = ftype;
204:
205: PetscFree((*B)->solvertype);
206: PetscStrallocpy(MATSOLVERBAS,&(*B)->solvertype);
207: return(0);
208: }