xref: /petsc/src/ksp/pc/impls/factor/qr/qr.c (revision a2fc1e05f8e2357324c93b1aa04eb0258ff33a90)
1*a2fc1e05SToby Isaac 
2*a2fc1e05SToby Isaac /*
3*a2fc1e05SToby Isaac    Defines a direct QR factorization preconditioner for any Mat implementation
4*a2fc1e05SToby Isaac    Note: this need not be considered a preconditioner since it supplies
5*a2fc1e05SToby Isaac          a direct solver.
6*a2fc1e05SToby Isaac */
7*a2fc1e05SToby Isaac #include <../src/ksp/pc/impls/factor/qr/qr.h>  /*I "petscpc.h" I*/
8*a2fc1e05SToby Isaac 
9*a2fc1e05SToby Isaac static PetscErrorCode PCSetUp_QR(PC pc)
10*a2fc1e05SToby Isaac {
11*a2fc1e05SToby Isaac   PetscErrorCode         ierr;
12*a2fc1e05SToby Isaac   PC_QR                  *dir = (PC_QR*)pc->data;
13*a2fc1e05SToby Isaac   MatSolverType          stype;
14*a2fc1e05SToby Isaac   MatFactorError         err;
15*a2fc1e05SToby Isaac 
16*a2fc1e05SToby Isaac   PetscFunctionBegin;
17*a2fc1e05SToby Isaac   pc->failedreason = PC_NOERROR;
18*a2fc1e05SToby Isaac   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
19*a2fc1e05SToby Isaac 
20*a2fc1e05SToby Isaac   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
21*a2fc1e05SToby Isaac   if (dir->hdr.inplace) {
22*a2fc1e05SToby Isaac     MatFactorType ftype;
23*a2fc1e05SToby Isaac 
24*a2fc1e05SToby Isaac     ierr = MatGetFactorType(pc->pmat, &ftype);CHKERRQ(ierr);
25*a2fc1e05SToby Isaac     if (ftype == MAT_FACTOR_NONE) {
26*a2fc1e05SToby Isaac       ierr = MatQRFactor(pc->pmat,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
27*a2fc1e05SToby Isaac       ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr);
28*a2fc1e05SToby Isaac       if (err) { /* Factor() fails */
29*a2fc1e05SToby Isaac         pc->failedreason = (PCFailedReason)err;
30*a2fc1e05SToby Isaac         PetscFunctionReturn(0);
31*a2fc1e05SToby Isaac       }
32*a2fc1e05SToby Isaac     }
33*a2fc1e05SToby Isaac     ((PC_Factor*)dir)->fact = pc->pmat;
34*a2fc1e05SToby Isaac   } else {
35*a2fc1e05SToby Isaac     MatInfo info;
36*a2fc1e05SToby Isaac 
37*a2fc1e05SToby Isaac     if (!pc->setupcalled) {
38*a2fc1e05SToby Isaac       if (!((PC_Factor*)dir)->fact) {
39*a2fc1e05SToby Isaac         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_QR,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
40*a2fc1e05SToby Isaac         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
41*a2fc1e05SToby Isaac       }
42*a2fc1e05SToby Isaac       ierr = MatQRFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
43*a2fc1e05SToby Isaac       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
44*a2fc1e05SToby Isaac       dir->hdr.actualfill = info.fill_ratio_needed;
45*a2fc1e05SToby Isaac     } else if (pc->flag != SAME_NONZERO_PATTERN) {
46*a2fc1e05SToby Isaac       ierr = MatQRFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
47*a2fc1e05SToby Isaac       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
48*a2fc1e05SToby Isaac       dir->hdr.actualfill = info.fill_ratio_needed;
49*a2fc1e05SToby Isaac     } else {
50*a2fc1e05SToby Isaac       ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
51*a2fc1e05SToby Isaac     }
52*a2fc1e05SToby Isaac     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
53*a2fc1e05SToby Isaac     if (err) { /* FactorSymbolic() fails */
54*a2fc1e05SToby Isaac       pc->failedreason = (PCFailedReason)err;
55*a2fc1e05SToby Isaac       PetscFunctionReturn(0);
56*a2fc1e05SToby Isaac     }
57*a2fc1e05SToby Isaac 
58*a2fc1e05SToby Isaac     ierr = MatQRFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
59*a2fc1e05SToby Isaac     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
60*a2fc1e05SToby Isaac     if (err) { /* FactorNumeric() fails */
61*a2fc1e05SToby Isaac       pc->failedreason = (PCFailedReason)err;
62*a2fc1e05SToby Isaac     }
63*a2fc1e05SToby Isaac   }
64*a2fc1e05SToby Isaac 
65*a2fc1e05SToby Isaac   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
66*a2fc1e05SToby Isaac   if (!stype) {
67*a2fc1e05SToby Isaac     MatSolverType solverpackage;
68*a2fc1e05SToby Isaac     ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr);
69*a2fc1e05SToby Isaac     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
70*a2fc1e05SToby Isaac   }
71*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
72*a2fc1e05SToby Isaac }
73*a2fc1e05SToby Isaac 
74*a2fc1e05SToby Isaac static PetscErrorCode PCReset_QR(PC pc)
75*a2fc1e05SToby Isaac {
76*a2fc1e05SToby Isaac   PC_QR          *dir = (PC_QR*)pc->data;
77*a2fc1e05SToby Isaac   PetscErrorCode ierr;
78*a2fc1e05SToby Isaac 
79*a2fc1e05SToby Isaac   PetscFunctionBegin;
80*a2fc1e05SToby Isaac   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
81*a2fc1e05SToby Isaac   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
82*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
83*a2fc1e05SToby Isaac }
84*a2fc1e05SToby Isaac 
85*a2fc1e05SToby Isaac static PetscErrorCode PCDestroy_QR(PC pc)
86*a2fc1e05SToby Isaac {
87*a2fc1e05SToby Isaac   PC_QR          *dir = (PC_QR*)pc->data;
88*a2fc1e05SToby Isaac   PetscErrorCode ierr;
89*a2fc1e05SToby Isaac 
90*a2fc1e05SToby Isaac   PetscFunctionBegin;
91*a2fc1e05SToby Isaac   ierr = PCReset_QR(pc);CHKERRQ(ierr);
92*a2fc1e05SToby Isaac   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
93*a2fc1e05SToby Isaac   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
94*a2fc1e05SToby Isaac   ierr = PetscFree(pc->data);CHKERRQ(ierr);
95*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
96*a2fc1e05SToby Isaac }
97*a2fc1e05SToby Isaac 
98*a2fc1e05SToby Isaac static PetscErrorCode PCApply_QR(PC pc,Vec x,Vec y)
99*a2fc1e05SToby Isaac {
100*a2fc1e05SToby Isaac   PC_QR          *dir = (PC_QR*)pc->data;
101*a2fc1e05SToby Isaac   Mat            fact;
102*a2fc1e05SToby Isaac   PetscErrorCode ierr;
103*a2fc1e05SToby Isaac 
104*a2fc1e05SToby Isaac   PetscFunctionBegin;
105*a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor*)dir)->fact;
106*a2fc1e05SToby Isaac   ierr = MatSolve(fact,x,y);CHKERRQ(ierr);
107*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
108*a2fc1e05SToby Isaac }
109*a2fc1e05SToby Isaac 
110*a2fc1e05SToby Isaac static PetscErrorCode PCMatApply_QR(PC pc,Mat X,Mat Y)
111*a2fc1e05SToby Isaac {
112*a2fc1e05SToby Isaac   PC_QR          *dir = (PC_QR*)pc->data;
113*a2fc1e05SToby Isaac   Mat            fact;
114*a2fc1e05SToby Isaac   PetscErrorCode ierr;
115*a2fc1e05SToby Isaac 
116*a2fc1e05SToby Isaac   PetscFunctionBegin;
117*a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor*)dir)->fact;
118*a2fc1e05SToby Isaac   ierr = MatMatSolve(fact,X,Y);CHKERRQ(ierr);
119*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
120*a2fc1e05SToby Isaac }
121*a2fc1e05SToby Isaac 
122*a2fc1e05SToby Isaac static PetscErrorCode PCApplyTranspose_QR(PC pc,Vec x,Vec y)
123*a2fc1e05SToby Isaac {
124*a2fc1e05SToby Isaac   PC_QR          *dir = (PC_QR*)pc->data;
125*a2fc1e05SToby Isaac   Mat            fact;
126*a2fc1e05SToby Isaac   PetscErrorCode ierr;
127*a2fc1e05SToby Isaac 
128*a2fc1e05SToby Isaac   PetscFunctionBegin;
129*a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor*)dir)->fact;
130*a2fc1e05SToby Isaac   ierr = MatSolveTranspose(fact,x,y);CHKERRQ(ierr);
131*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
132*a2fc1e05SToby Isaac }
133*a2fc1e05SToby Isaac 
134*a2fc1e05SToby Isaac /* -----------------------------------------------------------------------------------*/
135*a2fc1e05SToby Isaac 
136*a2fc1e05SToby Isaac /*MC
137*a2fc1e05SToby Isaac    PCQR - Uses a direct solver, based on QR factorization, as a preconditioner
138*a2fc1e05SToby Isaac 
139*a2fc1e05SToby Isaac    Level: beginner
140*a2fc1e05SToby Isaac 
141*a2fc1e05SToby Isaac    Notes:
142*a2fc1e05SToby Isaac     Usually this will compute an "exact" solution in one iteration and does
143*a2fc1e05SToby Isaac           not need a Krylov method (i.e. you can use -ksp_type preonly, or
144*a2fc1e05SToby Isaac           KSPSetType(ksp,KSPPREONLY) for the Krylov method
145*a2fc1e05SToby Isaac 
146*a2fc1e05SToby Isaac .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
147*a2fc1e05SToby Isaac            PCILU, PCLU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
148*a2fc1e05SToby Isaac            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
149*a2fc1e05SToby Isaac            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
150*a2fc1e05SToby Isaac            PCFactorReorderForNonzeroDiagonal()
151*a2fc1e05SToby Isaac M*/
152*a2fc1e05SToby Isaac 
153*a2fc1e05SToby Isaac PETSC_EXTERN PetscErrorCode PCCreate_QR(PC pc)
154*a2fc1e05SToby Isaac {
155*a2fc1e05SToby Isaac   PetscErrorCode ierr;
156*a2fc1e05SToby Isaac   PC_QR          *dir;
157*a2fc1e05SToby Isaac 
158*a2fc1e05SToby Isaac   PetscFunctionBegin;
159*a2fc1e05SToby Isaac   ierr     = PetscNewLog(pc,&dir);CHKERRQ(ierr);
160*a2fc1e05SToby Isaac   pc->data = (void*)dir;
161*a2fc1e05SToby Isaac   ierr     = PCFactorInitialize(pc, MAT_FACTOR_QR);CHKERRQ(ierr);
162*a2fc1e05SToby Isaac 
163*a2fc1e05SToby Isaac   dir->col                   = NULL;
164*a2fc1e05SToby Isaac   pc->ops->reset             = PCReset_QR;
165*a2fc1e05SToby Isaac   pc->ops->destroy           = PCDestroy_QR;
166*a2fc1e05SToby Isaac   pc->ops->apply             = PCApply_QR;
167*a2fc1e05SToby Isaac   pc->ops->matapply          = PCMatApply_QR;
168*a2fc1e05SToby Isaac   pc->ops->applytranspose    = PCApplyTranspose_QR;
169*a2fc1e05SToby Isaac   pc->ops->setup             = PCSetUp_QR;
170*a2fc1e05SToby Isaac   pc->ops->setfromoptions    = PCSetFromOptions_Factor;
171*a2fc1e05SToby Isaac   pc->ops->view              = PCView_Factor;
172*a2fc1e05SToby Isaac   pc->ops->applyrichardson   = NULL;
173*a2fc1e05SToby Isaac   PetscFunctionReturn(0);
174*a2fc1e05SToby Isaac }
175