xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 55ba2a515af5dbbcf670f4e560eb9724ec139bf1)
1 #define PETSCKSP_DLL
2 
3 /*
4    Defines a direct factorization preconditioner for any Mat implementation
5    Note: this need not be consided a preconditioner since it supplies
6          a direct solver.
7 */
8 #include "private/pcimpl.h"                /*I "petscpc.h" I*/
9 
10 typedef struct {
11   Mat             fact;             /* factored matrix */
12   PetscReal       actualfill;       /* actual fill in factor */
13   PetscTruth      inplace;          /* flag indicating in-place factorization */
14   IS              row,col;          /* index sets used for reordering */
15   MatOrderingType ordering;         /* matrix ordering */
16   PetscTruth      reuseordering;    /* reuses previous reordering computed */
17   PetscTruth      reusefill;        /* reuse fill from previous Cholesky */
18   MatFactorInfo   info;
19 } PC_Cholesky;
20 
21 EXTERN_C_BEGIN
22 #undef __FUNCT__
23 #define __FUNCT__ "PCFactorSetZeroPivot_Cholesky"
24 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Cholesky(PC pc,PetscReal z)
25 {
26   PC_Cholesky *ch;
27 
28   PetscFunctionBegin;
29   ch                 = (PC_Cholesky*)pc->data;
30   ch->info.zeropivot = z;
31   PetscFunctionReturn(0);
32 }
33 EXTERN_C_END
34 
35 EXTERN_C_BEGIN
36 #undef __FUNCT__
37 #define __FUNCT__ "PCFactorSetShiftNonzero_Cholesky"
38 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Cholesky(PC pc,PetscReal shift)
39 {
40   PC_Cholesky *dir;
41 
42   PetscFunctionBegin;
43   dir = (PC_Cholesky*)pc->data;
44   if (shift == (PetscReal) PETSC_DECIDE) {
45     dir->info.shiftnz = 1.e-12;
46   } else {
47     dir->info.shiftnz = shift;
48   }
49   PetscFunctionReturn(0);
50 }
51 EXTERN_C_END
52 
53 EXTERN_C_BEGIN
54 #undef __FUNCT__
55 #define __FUNCT__ "PCFactorSetShiftPd_Cholesky"
56 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Cholesky(PC pc,PetscTruth shift)
57 {
58   PC_Cholesky *dir;
59 
60   PetscFunctionBegin;
61   dir = (PC_Cholesky*)pc->data;
62   if (shift) {
63     dir->info.shift_fraction = 0.0;
64     dir->info.shiftpd = 1.0;
65   } else {
66     dir->info.shiftpd = 0.0;
67   }
68   PetscFunctionReturn(0);
69 }
70 EXTERN_C_END
71 
72 EXTERN_C_BEGIN
73 #undef __FUNCT__
74 #define __FUNCT__ "PCCholeskySetReuseOrdering_Cholesky"
75 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
76 {
77   PC_Cholesky *lu;
78 
79   PetscFunctionBegin;
80   lu               = (PC_Cholesky*)pc->data;
81   lu->reuseordering = flag;
82   PetscFunctionReturn(0);
83 }
84 EXTERN_C_END
85 
86 EXTERN_C_BEGIN
87 #undef __FUNCT__
88 #define __FUNCT__ "PCCholeskySetReuseFill_Cholesky"
89 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetReuseFill_Cholesky(PC pc,PetscTruth flag)
90 {
91   PC_Cholesky *lu;
92 
93   PetscFunctionBegin;
94   lu = (PC_Cholesky*)pc->data;
95   lu->reusefill = flag;
96   PetscFunctionReturn(0);
97 }
98 EXTERN_C_END
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "PCSetFromOptions_Cholesky"
102 static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
103 {
104   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
105   PetscErrorCode ierr;
106   PetscTruth     flg;
107   char           tname[256];
108   PetscFList     ordlist;
109 
110   PetscFunctionBegin;
111   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
112   ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr);
113   ierr = PetscOptionsName("-pc_cholesky_in_place","Form Cholesky in the same memory as the matrix","PCCholeskySetUseInPlace",&flg);CHKERRQ(ierr);
114   if (flg) {
115     ierr = PCCholeskySetUseInPlace(pc);CHKERRQ(ierr);
116   }
117   ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr);
118 
119   ierr = PetscOptionsName("-pc_cholesky_reuse_fill","Use fill from previous factorization","PCCholeskySetReuseFill",&flg);CHKERRQ(ierr);
120   if (flg) {
121     ierr = PCCholeskySetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
122   }
123   ierr = PetscOptionsName("-pc_cholesky_reuse_ordering","Reuse ordering from previous factorization","PCCholeskySetReuseOrdering",&flg);CHKERRQ(ierr);
124   if (flg) {
125     ierr = PCCholeskySetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
126   }
127 
128   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
129   ierr = PetscOptionsList("-pc_cholesky_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCCholeskySetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
130   if (flg) {
131     ierr = PCCholeskySetMatOrdering(pc,tname);CHKERRQ(ierr);
132   }
133   ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
134   if (flg) {
135     ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
136   }
137   ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
138   ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
139   if (flg) {
140     ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
141   }
142   ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
143 
144   ierr = PetscOptionsTail();CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "PCView_Cholesky"
150 static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
151 {
152   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
153   PetscErrorCode ierr;
154   PetscTruth     iascii,isstring;
155 
156   PetscFunctionBegin;
157   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
158   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
159   if (iascii) {
160     MatInfo info;
161 
162     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");CHKERRQ(ierr);}
163     else             {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");CHKERRQ(ierr);}
164     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
165     if (lu->fact) {
166       ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
167       ierr = PetscViewerASCIIPrintf(viewer,"    Cholesky nonzeros %G\n",info.nz_used);CHKERRQ(ierr);
168       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr);
169       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
170       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
171     }
172     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
173     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
174   } else if (isstring) {
175     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
176   } else {
177     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCholesky",((PetscObject)viewer)->type_name);
178   }
179   PetscFunctionReturn(0);
180 }
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "PCGetFactoredMatrix_Cholesky"
184 static PetscErrorCode PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
185 {
186   PC_Cholesky *dir = (PC_Cholesky*)pc->data;
187 
188   PetscFunctionBegin;
189   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
190   *mat = dir->fact;
191   PetscFunctionReturn(0);
192 }
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "PCSetUp_Cholesky"
196 static PetscErrorCode PCSetUp_Cholesky(PC pc)
197 {
198   PetscErrorCode ierr;
199   PetscTruth     flg;
200   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
201 
202   PetscFunctionBegin;
203   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
204 
205   if (dir->inplace) {
206     if (dir->row && dir->col && (dir->row != dir->col)) {
207       ierr = ISDestroy(dir->row);CHKERRQ(ierr);
208       dir->row = 0;
209     }
210     if (dir->col) {
211       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
212       dir->col = 0;
213     }
214     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
215     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
216       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
217       dir->col=0;
218     }
219     if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
220     ierr = MatCholeskyFactor(pc->pmat,dir->row,&dir->info);CHKERRQ(ierr);
221     dir->fact = pc->pmat;
222   } else {
223     MatInfo info;
224     if (!pc->setupcalled) {
225       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
226       if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
227         ierr = ISDestroy(dir->col);CHKERRQ(ierr);
228         dir->col=0;
229       }
230       ierr = PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
231       if (flg) {
232         PetscReal tol = 1.e-10;
233         ierr = PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
234         ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
235       }
236       if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
237       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
238       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
239       dir->actualfill = info.fill_ratio_needed;
240       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
241     } else if (pc->flag != SAME_NONZERO_PATTERN) {
242       if (!dir->reuseordering) {
243         if (dir->row && dir->col && (dir->row != dir->col)) {
244           ierr = ISDestroy(dir->row);CHKERRQ(ierr);
245           dir->row = 0;
246         }
247         if (dir->col) {
248           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
249           dir->col =0;
250         }
251         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
252         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
253           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
254           dir->col=0;
255         }
256         ierr = PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
257         if (flg) {
258           PetscReal tol = 1.e-10;
259           ierr = PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
260           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
261         }
262         if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
263       }
264       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
265       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
266       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
267       dir->actualfill = info.fill_ratio_needed;
268       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
269     }
270     ierr = MatCholeskyFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
271   }
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNCT__
276 #define __FUNCT__ "PCDestroy_Cholesky"
277 static PetscErrorCode PCDestroy_Cholesky(PC pc)
278 {
279   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
280   PetscErrorCode ierr;
281 
282   PetscFunctionBegin;
283   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
284   if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
285   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
286   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
287   ierr = PetscFree(dir);CHKERRQ(ierr);
288   PetscFunctionReturn(0);
289 }
290 
291 #undef __FUNCT__
292 #define __FUNCT__ "PCApply_Cholesky"
293 static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
294 {
295   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
296   PetscErrorCode ierr;
297 
298   PetscFunctionBegin;
299   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
300   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "PCApplyTranspose_Cholesky"
306 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
307 {
308   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
309   PetscErrorCode ierr;
310 
311   PetscFunctionBegin;
312   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
313   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
314   PetscFunctionReturn(0);
315 }
316 
317 /* -----------------------------------------------------------------------------------*/
318 
319 EXTERN_C_BEGIN
320 #undef __FUNCT__
321 #define __FUNCT__ "PCFactorSetFill_Cholesky"
322 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Cholesky(PC pc,PetscReal fill)
323 {
324   PC_Cholesky *dir;
325 
326   PetscFunctionBegin;
327   dir = (PC_Cholesky*)pc->data;
328   dir->info.fill = fill;
329   PetscFunctionReturn(0);
330 }
331 EXTERN_C_END
332 
333 EXTERN_C_BEGIN
334 #undef __FUNCT__
335 #define __FUNCT__ "PCCholeskySetUseInPlace_Cholesky"
336 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetUseInPlace_Cholesky(PC pc)
337 {
338   PC_Cholesky *dir;
339 
340   PetscFunctionBegin;
341   dir = (PC_Cholesky*)pc->data;
342   dir->inplace = PETSC_TRUE;
343   PetscFunctionReturn(0);
344 }
345 EXTERN_C_END
346 
347 EXTERN_C_BEGIN
348 #undef __FUNCT__
349 #define __FUNCT__ "PCCholeskySetMatOrdering_Cholesky"
350 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
351 {
352   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
353   PetscErrorCode ierr;
354 
355   PetscFunctionBegin;
356   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
357   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
358   PetscFunctionReturn(0);
359 }
360 EXTERN_C_END
361 
362 /* -----------------------------------------------------------------------------------*/
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "PCCholeskySetReuseOrdering"
366 /*@
367    PCCholeskySetReuseOrdering - When similar matrices are factored, this
368    causes the ordering computed in the first factor to be used for all
369    following factors.
370 
371    Collective on PC
372 
373    Input Parameters:
374 +  pc - the preconditioner context
375 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
376 
377    Options Database Key:
378 .  -pc_cholesky_reuse_ordering - Activate PCCholeskySetReuseOrdering()
379 
380    Level: intermediate
381 
382 .keywords: PC, levels, reordering, factorization, incomplete, LU
383 
384 .seealso: PCCholeskySetReuseFill(), PCICholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
385 @*/
386 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetReuseOrdering(PC pc,PetscTruth flag)
387 {
388   PetscErrorCode ierr,(*f)(PC,PetscTruth);
389 
390   PetscFunctionBegin;
391   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
392   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
393   if (f) {
394     ierr = (*f)(pc,flag);CHKERRQ(ierr);
395   }
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "PCCholeskySetReuseFill"
401 /*@
402    PCCholeskySetReuseFill - When matrices with same nonzero structure are Cholesky factored,
403    this causes later ones to use the fill computed in the initial factorization.
404 
405    Collective on PC
406 
407    Input Parameters:
408 +  pc - the preconditioner context
409 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
410 
411    Options Database Key:
412 .  -pc_cholesky_reuse_fill - Activates PCCholeskySetReuseFill()
413 
414    Level: intermediate
415 
416 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
417 
418 .seealso: PCICholeskySetReuseOrdering(), PCCholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
419 @*/
420 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetReuseFill(PC pc,PetscTruth flag)
421 {
422   PetscErrorCode ierr,(*f)(PC,PetscTruth);
423 
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
426   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
427   if (f) {
428     ierr = (*f)(pc,flag);CHKERRQ(ierr);
429   }
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "PCCholeskySetUseInPlace"
435 /*@
436    PCCholeskySetUseInPlace - Tells the system to do an in-place factorization.
437    For dense matrices, this enables the solution of much larger problems.
438    For sparse matrices the factorization cannot be done truly in-place
439    so this does not save memory during the factorization, but after the matrix
440    is factored, the original unfactored matrix is freed, thus recovering that
441    space.
442 
443    Collective on PC
444 
445    Input Parameters:
446 .  pc - the preconditioner context
447 
448    Options Database Key:
449 .  -pc_cholesky_in_place - Activates in-place factorization
450 
451    Notes:
452    PCCholeskySetUseInplace() can only be used with the KSP method KSPPREONLY or when
453    a different matrix is provided for the multiply and the preconditioner in
454    a call to KSPSetOperators().
455    This is because the Krylov space methods require an application of the
456    matrix multiplication, which is not possible here because the matrix has
457    been factored in-place, replacing the original matrix.
458 
459    Level: intermediate
460 
461 .keywords: PC, set, factorization, direct, inplace, in-place, Cholesky
462 
463 .seealso: PCICholeskySetUseInPlace()
464 @*/
465 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetUseInPlace(PC pc)
466 {
467   PetscErrorCode ierr,(*f)(PC);
468 
469   PetscFunctionBegin;
470   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
471   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
472   if (f) {
473     ierr = (*f)(pc);CHKERRQ(ierr);
474   }
475   PetscFunctionReturn(0);
476 }
477 
478 #undef __FUNCT__
479 #define __FUNCT__ "PCCholeskySetMatOrdering"
480 /*@
481     PCCholeskySetMatOrdering - Sets the ordering routine (to reduce fill) to
482     be used it the Cholesky factorization.
483 
484     Collective on PC
485 
486     Input Parameters:
487 +   pc - the preconditioner context
488 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
489 
490     Options Database Key:
491 .   -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
492 
493     Level: intermediate
494 
495 .seealso: PCICholeskySetMatOrdering()
496 @*/
497 PetscErrorCode PETSCKSP_DLLEXPORT PCCholeskySetMatOrdering(PC pc,MatOrderingType ordering)
498 {
499   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
500 
501   PetscFunctionBegin;
502   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
503   if (f) {
504     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
505   }
506   PetscFunctionReturn(0);
507 }
508 
509 /*MC
510    PCCholesky - Uses a direct solver, based on Cholesky factorization, as a preconditioner
511 
512    Options Database Keys:
513 +  -pc_cholesky_reuse_ordering - Activate PCLUSetReuseOrdering()
514 .  -pc_cholesky_reuse_fill - Activates PCLUSetReuseFill()
515 .  -pc_factor_fill <fill> - Sets fill amount
516 .  -pc_cholesky_in_place - Activates in-place factorization
517 .  -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
518 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
519 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
520    is optional with PETSC_TRUE being the default
521 
522    Notes: Not all options work for all matrix formats
523 
524    Level: beginner
525 
526    Concepts: Cholesky factorization, direct solver
527 
528    Notes: Usually this will compute an "exact" solution in one iteration and does
529           not need a Krylov method (i.e. you can use -ksp_type preonly, or
530           KSPSetType(ksp,KSPPREONLY) for the Krylov method
531 
532 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
533            PCILU, PCLU, PCICC, PCCholeskySetReuseOrdering(), PCCholeskySetReuseFill(), PCGetFactoredMatrix(),
534            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
535 	   PCCholeskySetUseInPlace(), PCCholeskySetMatOrdering(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
536 
537 M*/
538 
539 EXTERN_C_BEGIN
540 #undef __FUNCT__
541 #define __FUNCT__ "PCCreate_Cholesky"
542 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc)
543 {
544   PetscErrorCode ierr;
545   PC_Cholesky    *dir;
546 
547   PetscFunctionBegin;
548   ierr = PetscNew(PC_Cholesky,&dir);CHKERRQ(ierr);
549   ierr = PetscLogObjectMemory(pc,sizeof(PC_Cholesky));CHKERRQ(ierr);
550 
551   dir->fact                   = 0;
552   dir->inplace                = PETSC_FALSE;
553   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
554   dir->info.fill              = 5.0;
555   dir->info.shiftnz           = 0.0;
556   dir->info.shiftpd           = 0.0; /* false */
557   dir->info.shift_fraction    = 0.0;
558   dir->info.pivotinblocks     = 1.0;
559   dir->col                    = 0;
560   dir->row                    = 0;
561   ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
562   dir->reusefill        = PETSC_FALSE;
563   dir->reuseordering    = PETSC_FALSE;
564   pc->data              = (void*)dir;
565 
566   pc->ops->destroy           = PCDestroy_Cholesky;
567   pc->ops->apply             = PCApply_Cholesky;
568   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
569   pc->ops->setup             = PCSetUp_Cholesky;
570   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
571   pc->ops->view              = PCView_Cholesky;
572   pc->ops->applyrichardson   = 0;
573   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_Cholesky;
574 
575   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Cholesky",
576                     PCFactorSetZeroPivot_Cholesky);CHKERRQ(ierr);
577   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Cholesky",
578                     PCFactorSetShiftNonzero_Cholesky);CHKERRQ(ierr);
579   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Cholesky",
580                     PCFactorSetShiftPd_Cholesky);CHKERRQ(ierr);
581 
582   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Cholesky",
583                     PCFactorSetFill_Cholesky);CHKERRQ(ierr);
584   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetUseInPlace_C","PCCholeskySetUseInPlace_Cholesky",
585                     PCCholeskySetUseInPlace_Cholesky);CHKERRQ(ierr);
586   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetMatOrdering_C","PCCholeskySetMatOrdering_Cholesky",
587                     PCCholeskySetMatOrdering_Cholesky);CHKERRQ(ierr);
588   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseOrdering_C","PCCholeskySetReuseOrdering_Cholesky",
589                     PCCholeskySetReuseOrdering_Cholesky);CHKERRQ(ierr);
590   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseFill_C","PCCholeskySetReuseFill_Cholesky",
591                     PCCholeskySetReuseFill_Cholesky);CHKERRQ(ierr);
592   PetscFunctionReturn(0);
593 }
594 EXTERN_C_END
595