14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 34b9ad928SBarry Smith 4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/ 51f6cc5b2SSatish Balay /*@C 654b2cd4bSJed Brown PCMGResidualDefault - Default routine to calculate the residual. 74b9ad928SBarry Smith 84b9ad928SBarry Smith Collective on Mat and Vec 94b9ad928SBarry Smith 104b9ad928SBarry Smith Input Parameters: 114b9ad928SBarry Smith + mat - the matrix 124b9ad928SBarry Smith . b - the right-hand-side 134b9ad928SBarry Smith - x - the approximate solution 144b9ad928SBarry Smith 154b9ad928SBarry Smith Output Parameter: 164b9ad928SBarry Smith . r - location to store the residual 174b9ad928SBarry Smith 18d0e4de75SBarry Smith Level: developer 194b9ad928SBarry Smith 2097177400SBarry Smith .seealso: PCMGSetResidual() 214b9ad928SBarry Smith @*/ 2254b2cd4bSJed Brown PetscErrorCode PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r) 234b9ad928SBarry Smith { 24dfbe8321SBarry Smith PetscErrorCode ierr; 254b9ad928SBarry Smith 264b9ad928SBarry Smith PetscFunctionBegin; 27f9426fe0SMark Adams ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr); 284b9ad928SBarry Smith PetscFunctionReturn(0); 294b9ad928SBarry Smith } 304b9ad928SBarry Smith 31f39d8e23SSatish Balay /*@ 3297177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 334b9ad928SBarry Smith 344b9ad928SBarry Smith Not Collective 354b9ad928SBarry Smith 364b9ad928SBarry Smith Input Parameter: 374b9ad928SBarry Smith . pc - the multigrid context 384b9ad928SBarry Smith 394b9ad928SBarry Smith Output Parameter: 404b9ad928SBarry Smith . ksp - the coarse grid solver context 414b9ad928SBarry Smith 424b9ad928SBarry Smith Level: advanced 434b9ad928SBarry Smith 4428668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother() 454b9ad928SBarry Smith @*/ 467087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 474b9ad928SBarry Smith { 48f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 49f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 504b9ad928SBarry Smith 514b9ad928SBarry Smith PetscFunctionBegin; 52c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 53f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 544b9ad928SBarry Smith PetscFunctionReturn(0); 554b9ad928SBarry Smith } 564b9ad928SBarry Smith 574b9ad928SBarry Smith /*@C 5897177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 594b9ad928SBarry Smith on the lth level. 604b9ad928SBarry Smith 61ad4df100SBarry Smith Logically Collective on PC and Mat 624b9ad928SBarry Smith 634b9ad928SBarry Smith Input Parameters: 644b9ad928SBarry Smith + pc - the multigrid context 654b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 66157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 67d0e4de75SBarry Smith previous one were provided then a default is used 684b9ad928SBarry Smith - mat - matrix associated with residual 694b9ad928SBarry Smith 704b9ad928SBarry Smith Level: advanced 714b9ad928SBarry Smith 7254b2cd4bSJed Brown .seealso: PCMGResidualDefault() 734b9ad928SBarry Smith @*/ 747087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 754b9ad928SBarry Smith { 76f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 77f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 78298cc208SBarry Smith PetscErrorCode ierr; 794b9ad928SBarry Smith 804b9ad928SBarry Smith PetscFunctionBegin; 81c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 82ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 832fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 8454b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 85f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 86298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 87f3fbd535SBarry Smith mglevels[l]->A = mat; 884b9ad928SBarry Smith PetscFunctionReturn(0); 894b9ad928SBarry Smith } 904b9ad928SBarry Smith 914b9ad928SBarry Smith /*@ 92aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 93bf5b2e24SBarry Smith interpolation from l-1 to the lth level 944b9ad928SBarry Smith 95ad4df100SBarry Smith Logically Collective on PC and Mat 964b9ad928SBarry Smith 974b9ad928SBarry Smith Input Parameters: 984b9ad928SBarry Smith + pc - the multigrid context 994b9ad928SBarry Smith . mat - the interpolation operator 100bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1014b9ad928SBarry Smith 1024b9ad928SBarry Smith Level: advanced 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith Notes: 1054b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1064b9ad928SBarry Smith for the same level. 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1094b9ad928SBarry Smith out from the matrix size which one it is. 1104b9ad928SBarry Smith 11197177400SBarry Smith .seealso: PCMGSetRestriction() 1124b9ad928SBarry Smith @*/ 1137087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1144b9ad928SBarry Smith { 115f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 116f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 117fccaa45eSBarry Smith PetscErrorCode ierr; 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith PetscFunctionBegin; 120c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 121ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 122ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 123c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1246bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 1252fa5cd67SKarl Rupp 126f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1274b9ad928SBarry Smith PetscFunctionReturn(0); 1284b9ad928SBarry Smith } 1294b9ad928SBarry Smith 1308a2c336bSFande Kong /*@ 131*95750439SFande Kong PCMGSetOperators - Sets operator and preconditioning matrix for lth level 1328a2c336bSFande Kong 1338a2c336bSFande Kong Logically Collective on PC and Mat 1348a2c336bSFande Kong 1358a2c336bSFande Kong Input Parameters: 1368a2c336bSFande Kong + pc - the multigrid context 1378a2c336bSFande Kong . Amat - the operator 1388a2c336bSFande Kong . pmat - the preconditioning operator 139*95750439SFande Kong - l - the level (0 is the coarsest) to supply 1408a2c336bSFande Kong 1418a2c336bSFande Kong Level: advanced 1428a2c336bSFande Kong 1438a2c336bSFande Kong .keywords: multigrid, set, interpolate, level 1448a2c336bSFande Kong 1458a2c336bSFande Kong .seealso: PCMGSetRestriction(), PCMGSetInterpolation() 1468a2c336bSFande Kong @*/ 1478a2c336bSFande Kong PetscErrorCode PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat) 148360ee056SFande Kong { 149360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 150360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 151360ee056SFande Kong PetscErrorCode ierr; 152360ee056SFande Kong 153360ee056SFande Kong PetscFunctionBegin; 154360ee056SFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1558a2c336bSFande Kong PetscValidHeaderSpecific(Amat,MAT_CLASSID,3); 1568a2c336bSFande Kong PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4); 157360ee056SFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1588a2c336bSFande Kong ierr = KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat);CHKERRQ(ierr); 159360ee056SFande Kong PetscFunctionReturn(0); 160360ee056SFande Kong } 161360ee056SFande Kong 162c88c5224SJed Brown /*@ 163c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 164c88c5224SJed Brown interpolation from l-1 to the lth level 165c88c5224SJed Brown 166c88c5224SJed Brown Logically Collective on PC 167c88c5224SJed Brown 168c88c5224SJed Brown Input Parameters: 169c88c5224SJed Brown + pc - the multigrid context 170c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 171c88c5224SJed Brown 172c88c5224SJed Brown Output Parameter: 1733ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 174c88c5224SJed Brown 175c88c5224SJed Brown Level: advanced 176c88c5224SJed Brown 177c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 178c88c5224SJed Brown @*/ 179c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 180c88c5224SJed Brown { 181c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 182c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 183c88c5224SJed Brown PetscErrorCode ierr; 184c88c5224SJed Brown 185c88c5224SJed Brown PetscFunctionBegin; 186c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1873ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 188ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 189ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 190c88c5224SJed Brown if (!mglevels[l]->interpolate) { 1915aa31b60SBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 192c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 193c88c5224SJed Brown } 1943ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 195c88c5224SJed Brown PetscFunctionReturn(0); 196c88c5224SJed Brown } 197c88c5224SJed Brown 1984b9ad928SBarry Smith /*@ 1998a2c336bSFande Kong PCMGGetInterpolations - Gets interpolation matrices for all levels (except level 0) 2008a2c336bSFande Kong 2018a2c336bSFande Kong Logically Collective on PC 2028a2c336bSFande Kong 2038a2c336bSFande Kong Input Parameters: 2048a2c336bSFande Kong + pc - the multigrid context 2058a2c336bSFande Kong 2068a2c336bSFande Kong Output Parameter: 2078a2c336bSFande Kong - num_levels - the number of levels 2088a2c336bSFande Kong . interpolations - the interpolation matrices (size of num_levels-1) 2098a2c336bSFande Kong 2108a2c336bSFande Kong Notes: 2118a2c336bSFande Kong No new matrices are created, and the interpolations are the references to the original ones 2128a2c336bSFande Kong 2138a2c336bSFande Kong Level: advanced 2148a2c336bSFande Kong 2158a2c336bSFande Kong .keywords: MG, get, multigrid, interpolation, level 2168a2c336bSFande Kong 2178a2c336bSFande Kong .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale(), PCMGGetInterpolation() 2188a2c336bSFande Kong @*/ 2198a2c336bSFande Kong PetscErrorCode PCMGGetInterpolations(PC pc,PetscInt *num_levels,Mat *interpolations[]) 2208a2c336bSFande Kong { 2218a2c336bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 2228a2c336bSFande Kong PC_MG_Levels **mglevels = mg->levels; 2238a2c336bSFande Kong Mat *mat; 2248a2c336bSFande Kong PetscInt l; 2258a2c336bSFande Kong PetscErrorCode ierr; 2268a2c336bSFande Kong 2278a2c336bSFande Kong PetscFunctionBegin; 2288a2c336bSFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2298a2c336bSFande Kong PetscValidPointer(num_levels,2); 2308a2c336bSFande Kong PetscValidPointer(interpolations,3); 2318a2c336bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 2328a2c336bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 2338a2c336bSFande Kong for (l=1; l< mg->nlevels; l++) { 2348a2c336bSFande Kong mat[l-1] = mglevels[l]->interpolate; 2358a2c336bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr); 2368a2c336bSFande Kong } 2378a2c336bSFande Kong *num_levels = mg->nlevels; 2388a2c336bSFande Kong *interpolations = mat; 2398a2c336bSFande Kong PetscFunctionReturn(0); 2408a2c336bSFande Kong } 2418a2c336bSFande Kong 2428a2c336bSFande Kong /*@ 2438a2c336bSFande Kong PCMGGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level) 2448a2c336bSFande Kong 2458a2c336bSFande Kong Logically Collective on PC 2468a2c336bSFande Kong 2478a2c336bSFande Kong Input Parameters: 2488a2c336bSFande Kong + pc - the multigrid context 2498a2c336bSFande Kong 2508a2c336bSFande Kong Output Parameter: 2518a2c336bSFande Kong - num_levels - the number of levels 2528a2c336bSFande Kong . coarseOperators - the coarse operator matrices (size of num_levels-1) 2538a2c336bSFande Kong 2548a2c336bSFande Kong Notes: 2558a2c336bSFande Kong No new matrices are created, and the coarse operator matrices are the references to the original ones 2568a2c336bSFande Kong 2578a2c336bSFande Kong Level: advanced 2588a2c336bSFande Kong 2598a2c336bSFande Kong .keywords: MG, get, multigrid, interpolation, level 2608a2c336bSFande Kong 2618a2c336bSFande Kong .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale(), PCMGGetInterpolation(), PCMGGetInterpolations() 2628a2c336bSFande Kong @*/ 2638a2c336bSFande Kong PetscErrorCode PCMGGetCoarseOperators(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 2648a2c336bSFande Kong { 2658a2c336bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 2668a2c336bSFande Kong PC_MG_Levels **mglevels = mg->levels; 2678a2c336bSFande Kong PetscInt l; 2688a2c336bSFande Kong Mat *mat; 2698a2c336bSFande Kong PetscErrorCode ierr; 2708a2c336bSFande Kong 2718a2c336bSFande Kong PetscFunctionBegin; 2728a2c336bSFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2738a2c336bSFande Kong PetscValidPointer(num_levels,2); 2748a2c336bSFande Kong PetscValidPointer(coarseOperators,3); 2758a2c336bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 2768a2c336bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 2778a2c336bSFande Kong for (l=0; l<mg->nlevels-1; l++) { 2788a2c336bSFande Kong ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr); 2798a2c336bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr); 2808a2c336bSFande Kong } 2818a2c336bSFande Kong *num_levels = mg->nlevels; 2828a2c336bSFande Kong *coarseOperators = mat; 2838a2c336bSFande Kong PetscFunctionReturn(0); 2848a2c336bSFande Kong } 2858a2c336bSFande Kong 2868a2c336bSFande Kong /*@ 287eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 2884b9ad928SBarry Smith from level l to l-1. 2894b9ad928SBarry Smith 290ad4df100SBarry Smith Logically Collective on PC and Mat 2914b9ad928SBarry Smith 2924b9ad928SBarry Smith Input Parameters: 2934b9ad928SBarry Smith + pc - the multigrid context 294c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 295c88c5224SJed Brown - mat - the restriction matrix 2964b9ad928SBarry Smith 2974b9ad928SBarry Smith Level: advanced 2984b9ad928SBarry Smith 2994b9ad928SBarry Smith Notes: 3004b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 3014b9ad928SBarry Smith for the same level. 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 3044b9ad928SBarry Smith out from the matrix size which one it is. 3054b9ad928SBarry Smith 306aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 307fccaa45eSBarry Smith is used. 308fccaa45eSBarry Smith 309eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGetSetInjection() 3104b9ad928SBarry Smith @*/ 3117087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 3124b9ad928SBarry Smith { 313fccaa45eSBarry Smith PetscErrorCode ierr; 314f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 315f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3164b9ad928SBarry Smith 3174b9ad928SBarry Smith PetscFunctionBegin; 318c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 319c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 320ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 321ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 322c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 3236bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 3242fa5cd67SKarl Rupp 325f3fbd535SBarry Smith mglevels[l]->restrct = mat; 3264b9ad928SBarry Smith PetscFunctionReturn(0); 3274b9ad928SBarry Smith } 3284b9ad928SBarry Smith 329c88c5224SJed Brown /*@ 330eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 331c88c5224SJed Brown from level l to l-1. 332c88c5224SJed Brown 333c88c5224SJed Brown Logically Collective on PC and Mat 334c88c5224SJed Brown 335c88c5224SJed Brown Input Parameters: 336c88c5224SJed Brown + pc - the multigrid context 337c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 338c88c5224SJed Brown 339c88c5224SJed Brown Output Parameter: 340c88c5224SJed Brown . mat - the restriction matrix 341c88c5224SJed Brown 342c88c5224SJed Brown Level: advanced 343c88c5224SJed Brown 344eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection() 345c88c5224SJed Brown @*/ 346c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 347c88c5224SJed Brown { 348c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 349c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 350c88c5224SJed Brown PetscErrorCode ierr; 351c88c5224SJed Brown 352c88c5224SJed Brown PetscFunctionBegin; 353c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3543ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 355ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 356ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 357c88c5224SJed Brown if (!mglevels[l]->restrct) { 358ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 359c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 360c88c5224SJed Brown } 3613ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 362c88c5224SJed Brown PetscFunctionReturn(0); 363c88c5224SJed Brown } 364c88c5224SJed Brown 36573250ac0SBarry Smith /*@ 36673250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 36773250ac0SBarry Smith 368c88c5224SJed Brown Logically Collective on PC and Vec 369c88c5224SJed Brown 370c88c5224SJed Brown Input Parameters: 371c88c5224SJed Brown + pc - the multigrid context 372c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 373c88c5224SJed Brown . rscale - the scaling 374c88c5224SJed Brown 375c88c5224SJed Brown Level: advanced 376c88c5224SJed Brown 377c88c5224SJed Brown Notes: 378eab52d2bSLawrence Mitchell When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. It is preferable to use PCMGSetInjection() to control moving primal vectors. 379c88c5224SJed Brown 380eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 381c88c5224SJed Brown @*/ 382c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 383c88c5224SJed Brown { 384c88c5224SJed Brown PetscErrorCode ierr; 385c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 386c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 387c88c5224SJed Brown 388c88c5224SJed Brown PetscFunctionBegin; 389c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 390ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 391ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 392c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 393c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 3942fa5cd67SKarl Rupp 395c88c5224SJed Brown mglevels[l]->rscale = rscale; 396c88c5224SJed Brown PetscFunctionReturn(0); 397c88c5224SJed Brown } 398c88c5224SJed Brown 399c88c5224SJed Brown /*@ 400c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 401c88c5224SJed Brown 402c88c5224SJed Brown Collective on PC 40373250ac0SBarry Smith 40473250ac0SBarry Smith Input Parameters: 40573250ac0SBarry Smith + pc - the multigrid context 40673250ac0SBarry Smith . rscale - the scaling 40773250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 40873250ac0SBarry Smith 40973250ac0SBarry Smith Level: advanced 41073250ac0SBarry Smith 41173250ac0SBarry Smith Notes: 412eab52d2bSLawrence Mitchell When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. It is preferable to use PCMGGetInjection() to control moving primal vectors. 41373250ac0SBarry Smith 414eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 41573250ac0SBarry Smith @*/ 416c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 41773250ac0SBarry Smith { 41873250ac0SBarry Smith PetscErrorCode ierr; 41973250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 42073250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 42173250ac0SBarry Smith 42273250ac0SBarry Smith PetscFunctionBegin; 42373250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 424ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 425ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 426c88c5224SJed Brown if (!mglevels[l]->rscale) { 427c88c5224SJed Brown Mat R; 428c88c5224SJed Brown Vec X,Y,coarse,fine; 429c88c5224SJed Brown PetscInt M,N; 430c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 4312a7a6963SBarry Smith ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 432c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 4332fa5cd67SKarl Rupp if (M < N) { 4342fa5cd67SKarl Rupp fine = X; 4352fa5cd67SKarl Rupp coarse = Y; 4362fa5cd67SKarl Rupp } else if (N < M) { 4372fa5cd67SKarl Rupp fine = Y; coarse = X; 438ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 439c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 440c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 441c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 442c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 443c88c5224SJed Brown mglevels[l]->rscale = coarse; 444c88c5224SJed Brown } 445c88c5224SJed Brown *rscale = mglevels[l]->rscale; 44673250ac0SBarry Smith PetscFunctionReturn(0); 44773250ac0SBarry Smith } 44873250ac0SBarry Smith 449f39d8e23SSatish Balay /*@ 450eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 451eab52d2bSLawrence Mitchell from level l to l-1. 452eab52d2bSLawrence Mitchell 453eab52d2bSLawrence Mitchell Logically Collective on PC and Mat 454eab52d2bSLawrence Mitchell 455eab52d2bSLawrence Mitchell Input Parameters: 456eab52d2bSLawrence Mitchell + pc - the multigrid context 457eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 458eab52d2bSLawrence Mitchell - mat - the injection matrix 459eab52d2bSLawrence Mitchell 460eab52d2bSLawrence Mitchell Level: advanced 461eab52d2bSLawrence Mitchell 462eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction() 463eab52d2bSLawrence Mitchell @*/ 464eab52d2bSLawrence Mitchell PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 465eab52d2bSLawrence Mitchell { 466eab52d2bSLawrence Mitchell PetscErrorCode ierr; 467eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 468eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 469eab52d2bSLawrence Mitchell 470eab52d2bSLawrence Mitchell PetscFunctionBegin; 471eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 472eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 473eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 474eab52d2bSLawrence Mitchell if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 475eab52d2bSLawrence Mitchell ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 476eab52d2bSLawrence Mitchell ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr); 477eab52d2bSLawrence Mitchell 478eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 479eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 480eab52d2bSLawrence Mitchell } 481eab52d2bSLawrence Mitchell 482eab52d2bSLawrence Mitchell /*@ 483eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 484eab52d2bSLawrence Mitchell from level l to l-1. 485eab52d2bSLawrence Mitchell 486eab52d2bSLawrence Mitchell Logically Collective on PC and Mat 487eab52d2bSLawrence Mitchell 488eab52d2bSLawrence Mitchell Input Parameters: 489eab52d2bSLawrence Mitchell + pc - the multigrid context 490eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 491eab52d2bSLawrence Mitchell 492eab52d2bSLawrence Mitchell Output Parameter: 49399a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 494eab52d2bSLawrence Mitchell 495eab52d2bSLawrence Mitchell Level: advanced 496eab52d2bSLawrence Mitchell 497eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction() 498eab52d2bSLawrence Mitchell @*/ 499eab52d2bSLawrence Mitchell PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 500eab52d2bSLawrence Mitchell { 501eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 502eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 503eab52d2bSLawrence Mitchell 504eab52d2bSLawrence Mitchell PetscFunctionBegin; 505eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 506eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat,3); 507eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 508eab52d2bSLawrence Mitchell if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 509eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 510eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 511eab52d2bSLawrence Mitchell } 512eab52d2bSLawrence Mitchell 513eab52d2bSLawrence Mitchell /*@ 51497177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 51597177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 51697177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 5174b9ad928SBarry Smith post-smoothing. 5184b9ad928SBarry Smith 5194b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5204b9ad928SBarry Smith 5214b9ad928SBarry Smith Input Parameters: 5224b9ad928SBarry Smith + pc - the multigrid context 5234b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5244b9ad928SBarry Smith 5254b9ad928SBarry Smith Ouput Parameters: 5264b9ad928SBarry Smith . ksp - the smoother 5274b9ad928SBarry Smith 52857420d5bSBarry Smith Notes: 52957420d5bSBarry Smith Once you have called this routine, you can call KSPSetOperators(ksp,...) on the resulting ksp to provide the operators for the smoother for this level. 53057420d5bSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 53157420d5bSBarry Smith and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 53257420d5bSBarry Smith 5334b9ad928SBarry Smith Level: advanced 5344b9ad928SBarry Smith 53528668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 5364b9ad928SBarry Smith @*/ 5377087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 5384b9ad928SBarry Smith { 539f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 540f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5414b9ad928SBarry Smith 5424b9ad928SBarry Smith PetscFunctionBegin; 543c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 544f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5454b9ad928SBarry Smith PetscFunctionReturn(0); 5464b9ad928SBarry Smith } 5474b9ad928SBarry Smith 548f39d8e23SSatish Balay /*@ 54997177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 5504b9ad928SBarry Smith coarse grid correction (post-smoother). 5514b9ad928SBarry Smith 5524b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5534b9ad928SBarry Smith 5544b9ad928SBarry Smith Input Parameters: 5554b9ad928SBarry Smith + pc - the multigrid context 5564b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5574b9ad928SBarry Smith 5584b9ad928SBarry Smith Ouput Parameters: 5594b9ad928SBarry Smith . ksp - the smoother 5604b9ad928SBarry Smith 5614b9ad928SBarry Smith Level: advanced 5624b9ad928SBarry Smith 56395452b02SPatrick Sanan Notes: 56495452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 56589cce641SBarry Smith set options on the pre smoother also 56689cce641SBarry Smith 56797177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 5684b9ad928SBarry Smith @*/ 5697087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 5704b9ad928SBarry Smith { 571f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 572f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 573dfbe8321SBarry Smith PetscErrorCode ierr; 574f69a0ea3SMatthew Knepley const char *prefix; 5754b9ad928SBarry Smith MPI_Comm comm; 5764b9ad928SBarry Smith 5774b9ad928SBarry Smith PetscFunctionBegin; 578c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5794b9ad928SBarry Smith /* 5804b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 5814b9ad928SBarry Smith Thus we check if a different one has already been allocated, 5824b9ad928SBarry Smith if not we allocate it. 5834b9ad928SBarry Smith */ 584ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 585f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 58619fd82e9SBarry Smith KSPType ksptype; 58719fd82e9SBarry Smith PCType pctype; 588336babb1SJed Brown PC ipc; 589336babb1SJed Brown PetscReal rtol,abstol,dtol; 590336babb1SJed Brown PetscInt maxits; 591336babb1SJed Brown KSPNormType normtype; 592f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 593f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 594336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 5953bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 596336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 597336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 598336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 599336babb1SJed Brown 600f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 601422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 602f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 603f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 604336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 605336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 606336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 6070059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 608336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 609336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 6103bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 611ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 6124b9ad928SBarry Smith } 613f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 6144b9ad928SBarry Smith PetscFunctionReturn(0); 6154b9ad928SBarry Smith } 6164b9ad928SBarry Smith 617f39d8e23SSatish Balay /*@ 61897177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 6194b9ad928SBarry Smith coarse grid correction (pre-smoother). 6204b9ad928SBarry Smith 6214b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 6224b9ad928SBarry Smith 6234b9ad928SBarry Smith Input Parameters: 6244b9ad928SBarry Smith + pc - the multigrid context 6254b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 6264b9ad928SBarry Smith 6274b9ad928SBarry Smith Ouput Parameters: 6284b9ad928SBarry Smith . ksp - the smoother 6294b9ad928SBarry Smith 6304b9ad928SBarry Smith Level: advanced 6314b9ad928SBarry Smith 63295452b02SPatrick Sanan Notes: 63395452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 63489cce641SBarry Smith set options on the post smoother also 63589cce641SBarry Smith 63697177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 6374b9ad928SBarry Smith @*/ 6387087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 6394b9ad928SBarry Smith { 640dfbe8321SBarry Smith PetscErrorCode ierr; 641f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 642f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6434b9ad928SBarry Smith 6444b9ad928SBarry Smith PetscFunctionBegin; 645c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 6464b9ad928SBarry Smith /* make sure smoother up and down are different */ 647c5eb9154SBarry Smith if (l) { 6480298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 649d8148a5aSMatthew Knepley } 650f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 6514b9ad928SBarry Smith PetscFunctionReturn(0); 6524b9ad928SBarry Smith } 6534b9ad928SBarry Smith 6544b9ad928SBarry Smith /*@ 655cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 6564b9ad928SBarry Smith 657ad4df100SBarry Smith Logically Collective on PC 6584b9ad928SBarry Smith 6594b9ad928SBarry Smith Input Parameters: 6604b9ad928SBarry Smith + pc - the multigrid context 661c1cbb1deSBarry Smith . l - the level (0 is coarsest) 662c1cbb1deSBarry Smith - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 6634b9ad928SBarry Smith 6644b9ad928SBarry Smith Level: advanced 6654b9ad928SBarry Smith 666c1cbb1deSBarry Smith .seealso: PCMGSetCycleType() 6674b9ad928SBarry Smith @*/ 668c1cbb1deSBarry Smith PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 6694b9ad928SBarry Smith { 670f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 671f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6724b9ad928SBarry Smith 6734b9ad928SBarry Smith PetscFunctionBegin; 674c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 675ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 676c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 677c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc,c,3); 678f3fbd535SBarry Smith mglevels[l]->cycles = c; 6794b9ad928SBarry Smith PetscFunctionReturn(0); 6804b9ad928SBarry Smith } 6814b9ad928SBarry Smith 6824b9ad928SBarry Smith /*@ 68397177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 684fccaa45eSBarry Smith on a particular level. 6854b9ad928SBarry Smith 686ad4df100SBarry Smith Logically Collective on PC and Vec 6874b9ad928SBarry Smith 6884b9ad928SBarry Smith Input Parameters: 6894b9ad928SBarry Smith + pc - the multigrid context 6904b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6914b9ad928SBarry Smith - c - the space 6924b9ad928SBarry Smith 6934b9ad928SBarry Smith Level: advanced 6944b9ad928SBarry Smith 69595452b02SPatrick Sanan Notes: 69695452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 697fccaa45eSBarry Smith 698fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 699fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 700fccaa45eSBarry Smith 70197177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 7024b9ad928SBarry Smith @*/ 7037087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 7044b9ad928SBarry Smith { 705fccaa45eSBarry Smith PetscErrorCode ierr; 706f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 707f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7084b9ad928SBarry Smith 7094b9ad928SBarry Smith PetscFunctionBegin; 710c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 711ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 712ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 713c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 7146bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 7152fa5cd67SKarl Rupp 716f3fbd535SBarry Smith mglevels[l]->b = c; 7174b9ad928SBarry Smith PetscFunctionReturn(0); 7184b9ad928SBarry Smith } 7194b9ad928SBarry Smith 7204b9ad928SBarry Smith /*@ 72197177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 722fccaa45eSBarry Smith particular level. 7234b9ad928SBarry Smith 724ad4df100SBarry Smith Logically Collective on PC and Vec 7254b9ad928SBarry Smith 7264b9ad928SBarry Smith Input Parameters: 7274b9ad928SBarry Smith + pc - the multigrid context 728251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 7294b9ad928SBarry Smith - c - the space 7304b9ad928SBarry Smith 7314b9ad928SBarry Smith Level: advanced 7324b9ad928SBarry Smith 73395452b02SPatrick Sanan Notes: 73495452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 735fccaa45eSBarry Smith 736fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 737fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 738fccaa45eSBarry Smith 73997177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 7404b9ad928SBarry Smith @*/ 7417087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 7424b9ad928SBarry Smith { 743fccaa45eSBarry Smith PetscErrorCode ierr; 744f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 745f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7464b9ad928SBarry Smith 7474b9ad928SBarry Smith PetscFunctionBegin; 748c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 749ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 750ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 751c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 7526bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 7532fa5cd67SKarl Rupp 754f3fbd535SBarry Smith mglevels[l]->x = c; 7554b9ad928SBarry Smith PetscFunctionReturn(0); 7564b9ad928SBarry Smith } 7574b9ad928SBarry Smith 7584b9ad928SBarry Smith /*@ 75997177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 760fccaa45eSBarry Smith particular level. 7614b9ad928SBarry Smith 762ad4df100SBarry Smith Logically Collective on PC and Vec 7634b9ad928SBarry Smith 7644b9ad928SBarry Smith Input Parameters: 7654b9ad928SBarry Smith + pc - the multigrid context 7664b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 7674b9ad928SBarry Smith - c - the space 7684b9ad928SBarry Smith 7694b9ad928SBarry Smith Level: advanced 7704b9ad928SBarry Smith 77195452b02SPatrick Sanan Notes: 77295452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 773fccaa45eSBarry Smith 774fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 775fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 776fccaa45eSBarry Smith 7774b9ad928SBarry Smith @*/ 7787087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 7794b9ad928SBarry Smith { 780fccaa45eSBarry Smith PetscErrorCode ierr; 781f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 782f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7834b9ad928SBarry Smith 7844b9ad928SBarry Smith PetscFunctionBegin; 785c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 786ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 787ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 788c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 7896bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 7902fa5cd67SKarl Rupp 791f3fbd535SBarry Smith mglevels[l]->r = c; 7924b9ad928SBarry Smith PetscFunctionReturn(0); 7934b9ad928SBarry Smith } 794