xref: /petsc/src/mat/impls/aij/seq/bas/basfactor.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1b2f1ab58SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h>
5b2f1ab58SBarry Smith 
6b2f1ab58SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
7b2f1ab58SBarry Smith {
8b2f1ab58SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
9b2f1ab58SBarry Smith   Mat_SeqSBAIJ   *b;
10ace3abfcSBarry Smith   PetscBool      perm_identity,missing;
11b2f1ab58SBarry Smith   PetscInt       reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
12b2f1ab58SBarry Smith   const PetscInt *rip,*riip;
13b2f1ab58SBarry Smith   PetscInt       j;
14b2f1ab58SBarry Smith   PetscInt       d;
15b2f1ab58SBarry Smith   PetscInt       ncols,*cols,*uj;
16b2f1ab58SBarry Smith   PetscReal      fill=info->fill,levels=info->levels;
17b2f1ab58SBarry Smith   IS             iperm;
18b2f1ab58SBarry Smith   spbas_matrix   Pattern_0, Pattern_P;
19b2f1ab58SBarry Smith 
20b2f1ab58SBarry Smith   PetscFunctionBegin;
21*5f80ce2aSJacob Faibussowitsch   PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %" PetscInt_FMT " columns %" PetscInt_FMT,A->rmap->n,A->cmap->n);
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMissingDiagonal(A,&missing,&d));
23*5f80ce2aSJacob Faibussowitsch   PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,d);
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISIdentity(perm,&perm_identity));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISInvertPermutation(perm,PETSC_DECIDE,&iperm));
26b2f1ab58SBarry Smith 
27b2f1ab58SBarry Smith   /* ICC(0) without matrix ordering: simply copies fill pattern */
28b2f1ab58SBarry Smith   if (!levels && perm_identity) {
29*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(am+1,&ui));
30b2f1ab58SBarry Smith     ui[0] = 0;
31b2f1ab58SBarry Smith 
32b2f1ab58SBarry Smith     for (i=0; i<am; i++) {
33b2f1ab58SBarry Smith       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
34b2f1ab58SBarry Smith     }
35*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(ui[am]+1,&uj));
36b2f1ab58SBarry Smith     cols = uj;
37b2f1ab58SBarry Smith     for (i=0; i<am; i++) {
38b2f1ab58SBarry Smith       aj    = a->j + a->diag[i];
39b2f1ab58SBarry Smith       ncols = ui[i+1] - ui[i];
40b2f1ab58SBarry Smith       for (j=0; j<ncols; j++) *cols++ = *aj++;
41b2f1ab58SBarry Smith     }
42b2f1ab58SBarry Smith   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(iperm,&riip));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(perm,&rip));
45b2f1ab58SBarry Smith 
4671d55d18SBarry Smith     /* Create spbas_matrix for pattern */
47*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_pattern_only(am, am, ai, aj, &Pattern_0));
48b2f1ab58SBarry Smith 
4971d55d18SBarry Smith     /* Apply the permutation */
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_apply_reordering(&Pattern_0, rip, riip));
51b2f1ab58SBarry Smith 
5271d55d18SBarry Smith     /* Raise the power */
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_power(Pattern_0, (int) levels+1, &Pattern_P));
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_delete(Pattern_0));
55b2f1ab58SBarry Smith 
5671d55d18SBarry Smith     /* Keep only upper triangle of pattern */
57*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_keep_upper(&Pattern_P));
58b2f1ab58SBarry Smith 
5971d55d18SBarry Smith     /* Convert to Sparse Row Storage  */
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_matrix_to_crs(Pattern_P, NULL, &ui, &uj));
61*5f80ce2aSJacob Faibussowitsch     CHKERRQ(spbas_delete(Pattern_P));
62b2f1ab58SBarry Smith   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
63b2f1ab58SBarry Smith 
64b2f1ab58SBarry Smith   /* put together the new matrix in MATSEQSBAIJ format */
65b2f1ab58SBarry Smith 
66b2f1ab58SBarry Smith   b               = (Mat_SeqSBAIJ*)(fact)->data;
67b2f1ab58SBarry Smith   b->singlemalloc = PETSC_FALSE;
682205254eSKarl Rupp 
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ui[am]+1,&b->a));
702205254eSKarl Rupp 
71b2f1ab58SBarry Smith   b->j    = uj;
72b2f1ab58SBarry Smith   b->i    = ui;
73f4259b30SLisandro Dalcin   b->diag = NULL;
74f4259b30SLisandro Dalcin   b->ilen = NULL;
75f4259b30SLisandro Dalcin   b->imax = NULL;
76b2f1ab58SBarry Smith   b->row  = perm;
77b2f1ab58SBarry Smith   b->col  = perm;
782205254eSKarl Rupp 
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
812205254eSKarl Rupp 
82b2f1ab58SBarry Smith   b->icol          = iperm;
83b2f1ab58SBarry Smith   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(am+1,&b->solve_work));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogObjectMemory((PetscObject)(fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar))));
86b2f1ab58SBarry Smith   b->maxnz         = b->nz = ui[am];
87b2f1ab58SBarry Smith   b->free_a        = PETSC_TRUE;
88b2f1ab58SBarry Smith   b->free_ij       = PETSC_TRUE;
89b2f1ab58SBarry Smith 
90b2f1ab58SBarry Smith   (fact)->info.factor_mallocs   = reallocs;
91b2f1ab58SBarry Smith   (fact)->info.fill_ratio_given = fill;
92b2f1ab58SBarry Smith   if (ai[am] != 0) {
93b2f1ab58SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
94b2f1ab58SBarry Smith   } else {
95b2f1ab58SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
96b2f1ab58SBarry Smith   }
97b2f1ab58SBarry Smith   /*  (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */
98b2f1ab58SBarry Smith   PetscFunctionReturn(0);
99b2f1ab58SBarry Smith }
100b2f1ab58SBarry Smith 
101b2f1ab58SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo *info)
102b2f1ab58SBarry Smith {
103b2f1ab58SBarry Smith   Mat            C = B;
104b2f1ab58SBarry Smith   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
105b2f1ab58SBarry Smith   IS             ip=b->row,iip = b->icol;
106b2f1ab58SBarry Smith   const PetscInt *rip,*riip;
107b2f1ab58SBarry Smith   PetscInt       mbs=A->rmap->n,*bi=b->i,*bj=b->j;
108b2f1ab58SBarry Smith   MatScalar      *ba     = b->a;
109c88783f4SHong Zhang   PetscReal      shiftnz = info->shiftamount;
110d36aa34eSBarry Smith   PetscReal      droptol = -1;
111ace3abfcSBarry Smith   PetscBool      perm_identity;
112b2f1ab58SBarry Smith   spbas_matrix   Pattern, matrix_L,matrix_LT;
1139767453cSBarry Smith   PetscReal      mem_reduction;
114b2f1ab58SBarry Smith 
115b2f1ab58SBarry Smith   PetscFunctionBegin;
11671d55d18SBarry Smith   /* Reduce memory requirements:   erase values of B-matrix */
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ba));
11871d55d18SBarry Smith   /*   Compress (maximum) sparseness pattern of B-matrix */
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(bi));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(bj));
122b2f1ab58SBarry Smith 
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL,"    compression rate for spbas_compress_pattern %g \n",(double)mem_reduction));
124b2f1ab58SBarry Smith 
1254e1ff37aSBarry Smith   /* Make Cholesky decompositions with larger Manteuffel shifts until no more    negative diagonals are found. */
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(ip,&rip));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iip,&riip));
128b2f1ab58SBarry Smith 
129*5f80ce2aSJacob Faibussowitsch   if (info->usedt) droptol = info->dt;
130*5f80ce2aSJacob Faibussowitsch 
131*5f80ce2aSJacob Faibussowitsch   for (PetscErrorCode ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;) {
132cc442ecdSBarry Smith     PetscBool success;
133*5f80ce2aSJacob Faibussowitsch 
134*5f80ce2aSJacob Faibussowitsch     ierr = spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT,&success);
135cc442ecdSBarry Smith     if (!success) {
136b2f1ab58SBarry Smith       shiftnz *= 1.5;
1379767453cSBarry Smith       if (shiftnz < 1e-5) shiftnz=1e-5;
138*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n",(double)shiftnz));
139b2f1ab58SBarry Smith     }
140b2f1ab58SBarry Smith   }
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(spbas_delete(Pattern));
142b2f1ab58SBarry Smith 
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL,"    memory_usage for  spbas_incomplete_cholesky  %g bytes per row\n", (double)(PetscReal) (spbas_memory_requirement(matrix_LT)/ (PetscReal) mbs)));
144b2f1ab58SBarry Smith 
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(ip,&rip));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iip,&riip));
147b2f1ab58SBarry Smith 
14871d55d18SBarry Smith   /* Convert spbas_matrix to compressed row storage */
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(spbas_transpose(matrix_LT, &matrix_L));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(spbas_delete(matrix_LT));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj));
152b2f1ab58SBarry Smith   b->i =bi; b->j=bj; b->a=ba;
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(spbas_delete(matrix_L));
154b2f1ab58SBarry Smith 
15571d55d18SBarry Smith   /* Set the appropriate solution functions */
156*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISIdentity(ip,&perm_identity));
157b2f1ab58SBarry Smith   if (perm_identity) {
158b2f1ab58SBarry Smith     (B)->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
159b2f1ab58SBarry Smith     (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
160b2f1ab58SBarry Smith     (B)->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
161b2f1ab58SBarry Smith     (B)->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
162b2f1ab58SBarry Smith   } else {
163b2f1ab58SBarry Smith     (B)->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
164b2f1ab58SBarry Smith     (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
165b2f1ab58SBarry Smith     (B)->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
166b2f1ab58SBarry Smith     (B)->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
167b2f1ab58SBarry Smith   }
168b2f1ab58SBarry Smith 
169b2f1ab58SBarry Smith   C->assembled    = PETSC_TRUE;
170b2f1ab58SBarry Smith   C->preallocated = PETSC_TRUE;
1712205254eSKarl Rupp 
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(C->rmap->n));
173b2f1ab58SBarry Smith   PetscFunctionReturn(0);
174b2f1ab58SBarry Smith }
175b2f1ab58SBarry Smith 
176ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_bas(Mat A,MatSolverType *type)
1773dfa2556SBarry Smith {
1783dfa2556SBarry Smith   PetscFunctionBegin;
1793dfa2556SBarry Smith   *type = MATSOLVERBAS;
1803dfa2556SBarry Smith   PetscFunctionReturn(0);
1813dfa2556SBarry Smith }
1823dfa2556SBarry Smith 
183cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B)
184b2f1ab58SBarry Smith {
185b2f1ab58SBarry Smith   PetscInt       n = A->rmap->n;
186b2f1ab58SBarry Smith 
187b2f1ab58SBarry Smith   PetscFunctionBegin;
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),B));
189*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(*B,n,n,n,n));
190b2f1ab58SBarry Smith   if (ftype == MAT_FACTOR_ICC) {
191*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(*B,MATSEQSBAIJ));
192*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,NULL));
1932205254eSKarl Rupp 
194b2f1ab58SBarry Smith     (*B)->ops->iccfactorsymbolic     = MatICCFactorSymbolic_SeqAIJ_Bas;
195b2f1ab58SBarry Smith     (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
196*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectComposeFunction((PetscObject)*B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_bas));
197*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrallocpy(MATORDERINGND,(char**)&(*B)->preferredordering[MAT_FACTOR_LU]));
198*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrallocpy(MATORDERINGND,(char**)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
199e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
200d5f3da31SBarry Smith   (*B)->factortype = ftype;
20100c67f3bSHong Zhang 
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*B)->solvertype));
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(MATSOLVERBAS,&(*B)->solvertype));
204f73b0415SBarry Smith   (*B)->canuseordering = PETSC_TRUE;
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(MATORDERINGNATURAL,(char**)&(*B)->preferredordering[MAT_FACTOR_ICC]));
206b2f1ab58SBarry Smith   PetscFunctionReturn(0);
207b2f1ab58SBarry Smith }
208