xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
19b54502bSHong Zhang 
29b54502bSHong Zhang /*
39b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
49b54502bSHong Zhang    Note: this need not be consided a preconditioner since it supplies
59b54502bSHong Zhang          a direct solver.
69b54502bSHong Zhang */
7c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/
89b54502bSHong Zhang 
99b54502bSHong Zhang typedef struct {
10075768bcSBarry Smith   PC_Factor hdr;
119b54502bSHong Zhang   IS        row, col; /* index sets used for reordering */
129b54502bSHong Zhang } PC_Cholesky;
139b54502bSHong Zhang 
149371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Cholesky(PC pc, PetscOptionItems *PetscOptionsObject) {
159b54502bSHong Zhang   PetscFunctionBegin;
16d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Cholesky options");
17dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
18d0609cedSBarry Smith   PetscOptionsHeadEnd();
199b54502bSHong Zhang   PetscFunctionReturn(0);
209b54502bSHong Zhang }
219b54502bSHong Zhang 
229371c9d4SSatish Balay static PetscErrorCode PCSetUp_Cholesky(PC pc) {
23ace3abfcSBarry Smith   PetscBool      flg;
249b54502bSHong Zhang   PC_Cholesky   *dir = (PC_Cholesky *)pc->data;
25ea799195SBarry Smith   MatSolverType  stype;
2600e125f8SBarry Smith   MatFactorError err;
27f023e1d5SPierre Jolivet   const char    *prefix;
289b54502bSHong Zhang 
299b54502bSHong Zhang   PetscFunctionBegin;
30c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
313d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
329b54502bSHong Zhang 
3326cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
3426cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
3526cc229bSBarry Smith 
369566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
373d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
3848a46eb9SPierre Jolivet     if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row));
399566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&dir->col));
402c7c0729SBarry Smith     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
419566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
429566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
439b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
449566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
459b54502bSHong Zhang     }
469566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info));
479566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(pc->pmat, &err));
4800e125f8SBarry Smith     if (err) { /* Factor() fails */
4900e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
506baea169SHong Zhang       PetscFunctionReturn(0);
516baea169SHong Zhang     }
522fa5cd67SKarl Rupp 
53075768bcSBarry Smith     ((PC_Factor *)dir)->fact = pc->pmat;
549b54502bSHong Zhang   } else {
559b54502bSHong Zhang     MatInfo info;
5600e125f8SBarry Smith 
579b54502bSHong Zhang     if (!pc->setupcalled) {
58f73b0415SBarry Smith       PetscBool canuseordering;
59*4dfa11a4SJacob Faibussowitsch       if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); }
609566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
61f73b0415SBarry Smith       if (canuseordering) {
629566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
639566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
649bfd6278SHong Zhang         /* check if dir->row == dir->col */
654ac6704cSBarry Smith         if (dir->row) {
669566063dSJacob Faibussowitsch           PetscCall(ISEqual(dir->row, dir->col, &flg));
6728b400f6SJacob Faibussowitsch           PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal");
684ac6704cSBarry Smith         }
699566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */
709bfd6278SHong Zhang 
7190d69ab7SBarry Smith         flg = PETSC_FALSE;
729566063dSJacob Faibussowitsch         PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL));
739b54502bSHong Zhang         if (flg) {
749b54502bSHong Zhang           PetscReal tol = 1.e-10;
759566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
769566063dSJacob Faibussowitsch           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
779b54502bSHong Zhang         }
78a1f19f5aSHong Zhang       }
799566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
809566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
813d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
829b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
833d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
84f73b0415SBarry Smith         PetscBool canuseordering;
859566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
869566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact));
879566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
88f73b0415SBarry Smith         if (canuseordering) {
899566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->row));
909566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
919566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
929566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */
939d61c4f5SHong Zhang 
9490d69ab7SBarry Smith           flg = PETSC_FALSE;
959566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL));
969b54502bSHong Zhang           if (flg) {
979b54502bSHong Zhang             PetscReal tol = 1.e-10;
989566063dSJacob Faibussowitsch             PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
999566063dSJacob Faibussowitsch             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
1009b54502bSHong Zhang           }
1019b54502bSHong Zhang         }
1022c7c0729SBarry Smith       }
1039566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
1049566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1053d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
10604545d6dSBarry Smith     } else {
1079566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
108160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1099566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
110b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
11104545d6dSBarry Smith       }
1129b54502bSHong Zhang     }
1139566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
11400e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
11500e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1166baea169SHong Zhang       PetscFunctionReturn(0);
1176baea169SHong Zhang     }
1186baea169SHong Zhang 
1199566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1209566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
12100e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
12200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1236baea169SHong Zhang     }
1249b54502bSHong Zhang   }
12500c67f3bSHong Zhang 
1269566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
12700c67f3bSHong Zhang   if (!stype) {
128ea799195SBarry Smith     MatSolverType solverpackage;
1299566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1309566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
13100c67f3bSHong Zhang   }
1329b54502bSHong Zhang   PetscFunctionReturn(0);
1339b54502bSHong Zhang }
1349b54502bSHong Zhang 
1359371c9d4SSatish Balay static PetscErrorCode PCReset_Cholesky(PC pc) {
1369b54502bSHong Zhang   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1379b54502bSHong Zhang 
1389b54502bSHong Zhang   PetscFunctionBegin;
1399566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->row));
1419566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
14269d2c0f9SBarry Smith   PetscFunctionReturn(0);
14369d2c0f9SBarry Smith }
14469d2c0f9SBarry Smith 
1459371c9d4SSatish Balay static PetscErrorCode PCDestroy_Cholesky(PC pc) {
14669d2c0f9SBarry Smith   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
14769d2c0f9SBarry Smith 
14869d2c0f9SBarry Smith   PetscFunctionBegin;
1499566063dSJacob Faibussowitsch   PetscCall(PCReset_Cholesky(pc));
1509566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
1522e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1539566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1549b54502bSHong Zhang   PetscFunctionReturn(0);
1559b54502bSHong Zhang }
1569b54502bSHong Zhang 
1579371c9d4SSatish Balay static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y) {
1589b54502bSHong Zhang   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1599b54502bSHong Zhang 
1609b54502bSHong Zhang   PetscFunctionBegin;
1613d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1629566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat, x, y));
1632fa5cd67SKarl Rupp   } else {
1649566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1652fa5cd67SKarl Rupp   }
1669b54502bSHong Zhang   PetscFunctionReturn(0);
1679b54502bSHong Zhang }
1689b54502bSHong Zhang 
1699371c9d4SSatish Balay static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y) {
1707b6e2003SPierre Jolivet   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1717b6e2003SPierre Jolivet 
1727b6e2003SPierre Jolivet   PetscFunctionBegin;
1737b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1749566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat, X, Y));
1757b6e2003SPierre Jolivet   } else {
1769566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1777b6e2003SPierre Jolivet   }
1787b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1797b6e2003SPierre Jolivet }
1807b6e2003SPierre Jolivet 
1819371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y) {
1823d72fe1eSJed Brown   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1833d72fe1eSJed Brown 
1843d72fe1eSJed Brown   PetscFunctionBegin;
1853d72fe1eSJed Brown   if (dir->hdr.inplace) {
1869566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(pc->pmat, x, y));
1873d72fe1eSJed Brown   } else {
1889566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y));
1893d72fe1eSJed Brown   }
1903d72fe1eSJed Brown   PetscFunctionReturn(0);
1913d72fe1eSJed Brown }
1923d72fe1eSJed Brown 
1939371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y) {
1943d72fe1eSJed Brown   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1953d72fe1eSJed Brown 
1963d72fe1eSJed Brown   PetscFunctionBegin;
1973d72fe1eSJed Brown   if (dir->hdr.inplace) {
1989566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(pc->pmat, x, y));
1993d72fe1eSJed Brown   } else {
2009566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y));
2013d72fe1eSJed Brown   }
2023d72fe1eSJed Brown   PetscFunctionReturn(0);
2033d72fe1eSJed Brown }
2043d72fe1eSJed Brown 
2059371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y) {
2069b54502bSHong Zhang   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2079b54502bSHong Zhang 
2089b54502bSHong Zhang   PetscFunctionBegin;
2093d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2109566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat, x, y));
2112fa5cd67SKarl Rupp   } else {
2129566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
2132fa5cd67SKarl Rupp   }
2149b54502bSHong Zhang   PetscFunctionReturn(0);
2159b54502bSHong Zhang }
2169b54502bSHong Zhang 
2179b54502bSHong Zhang /*@
2182401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2199b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2209b54502bSHong Zhang    following factors.
2219b54502bSHong Zhang 
222f1580f4eSBarry Smith    Logically Collective on pc
2239b54502bSHong Zhang 
2249b54502bSHong Zhang    Input Parameters:
2259b54502bSHong Zhang +  pc - the preconditioner context
226f1580f4eSBarry Smith -  flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`
2279b54502bSHong Zhang 
2289b54502bSHong Zhang    Options Database Key:
229f1580f4eSBarry Smith .  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
2309b54502bSHong Zhang 
2319b54502bSHong Zhang    Level: intermediate
2329b54502bSHong Zhang 
233f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()`
2349b54502bSHong Zhang @*/
2359371c9d4SSatish Balay PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag) {
2369b54502bSHong Zhang   PetscFunctionBegin;
2370700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
238acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, flag, 2);
239cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag));
2409b54502bSHong Zhang   PetscFunctionReturn(0);
2419b54502bSHong Zhang }
2429b54502bSHong Zhang 
2439b54502bSHong Zhang /*MC
24496fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2459b54502bSHong Zhang 
2469b54502bSHong Zhang    Options Database Keys:
247f1580f4eSBarry Smith +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
248f1580f4eSBarry Smith .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
249f1580f4eSBarry Smith .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
25055ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2512401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
252145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2539b54502bSHong Zhang 
2549b54502bSHong Zhang    Level: beginner
2559b54502bSHong Zhang 
25695452b02SPatrick Sanan    Notes:
257f1580f4eSBarry Smith    Not all options work for all matrix formats
258f1580f4eSBarry Smith 
25995452b02SPatrick Sanan    Usually this will compute an "exact" solution in one iteration and does
2609b54502bSHong Zhang    not need a Krylov method (i.e. you can use -ksp_type preonly, or
261f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
2629b54502bSHong Zhang 
263db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
264db781477SPatrick Sanan           `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
265db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
266f1580f4eSBarry Smith           `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`
2679b54502bSHong Zhang M*/
2689b54502bSHong Zhang 
2699371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) {
2709b54502bSHong Zhang   PC_Cholesky *dir;
2719b54502bSHong Zhang 
2729b54502bSHong Zhang   PetscFunctionBegin;
273*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dir));
2743d1c1ea0SBarry Smith   pc->data = (void *)dir;
2759566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY));
2762fa5cd67SKarl Rupp 
277075768bcSBarry Smith   ((PC_Factor *)dir)->info.fill = 5.0;
2782fa5cd67SKarl Rupp 
2799b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
28069d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
2819b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
2827b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Cholesky;
2833d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
2843d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
2859b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
2869b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
2879b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
28892e08861SBarry Smith   pc->ops->view                = PCView_Factor;
2890a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
2909b54502bSHong Zhang   PetscFunctionReturn(0);
2919b54502bSHong Zhang }
292