xref: /petsc/src/ksp/pc/impls/kaczmarz/kaczmarz.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
23f93e5bdSPeter Brune 
33f93e5bdSPeter Brune typedef struct {
43f93e5bdSPeter Brune   PetscReal lambda;    /* damping parameter */
569b97631SPeter Brune   PetscBool symmetric; /* apply the projections symmetrically */
63f93e5bdSPeter Brune } PC_Kaczmarz;
73f93e5bdSPeter Brune 
PCDestroy_Kaczmarz(PC pc)8d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Kaczmarz(PC pc)
9d71ae5a4SJacob Faibussowitsch {
103f93e5bdSPeter Brune   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133f93e5bdSPeter Brune }
143f93e5bdSPeter Brune 
PCApply_Kaczmarz(PC pc,Vec x,Vec y)15d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Kaczmarz(PC pc, Vec x, Vec y)
16d71ae5a4SJacob Faibussowitsch {
173f93e5bdSPeter Brune   PC_Kaczmarz       *jac = (PC_Kaczmarz *)pc->data;
183f93e5bdSPeter Brune   PetscInt           xs, xe, ys, ye, ncols, i, j;
193f93e5bdSPeter Brune   const PetscInt    *cols;
20e298bee8SBarry Smith   const PetscScalar *vals, *xarray;
213f93e5bdSPeter Brune   PetscScalar        r;
223f93e5bdSPeter Brune   PetscReal          anrm;
23e298bee8SBarry Smith   PetscScalar       *yarray;
243f93e5bdSPeter Brune   PetscReal          lambda = jac->lambda;
253f93e5bdSPeter Brune 
263f93e5bdSPeter Brune   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(pc->pmat, &xs, &xe));
289566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(pc->pmat, &ys, &ye));
299566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.));
309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
323f93e5bdSPeter Brune   for (i = xs; i < xe; i++) {
333f93e5bdSPeter Brune     /* get the maximum row width and row norms */
349566063dSJacob Faibussowitsch     PetscCall(MatGetRow(pc->pmat, i, &ncols, &cols, &vals));
353f93e5bdSPeter Brune     r    = xarray[i - xs];
363f93e5bdSPeter Brune     anrm = 0.;
373f93e5bdSPeter Brune     for (j = 0; j < ncols; j++) {
38ad540459SPierre Jolivet       if (cols[j] >= ys && cols[j] < ye) r -= yarray[cols[j] - ys] * vals[j];
3969b97631SPeter Brune       anrm += PetscRealPart(PetscSqr(vals[j]));
403f93e5bdSPeter Brune     }
4169b97631SPeter Brune     if (anrm > 0.) {
423f93e5bdSPeter Brune       for (j = 0; j < ncols; j++) {
43ad540459SPierre Jolivet         if (cols[j] >= ys && cols[j] < ye) yarray[cols[j] - ys] += vals[j] * lambda * r / anrm;
443f93e5bdSPeter Brune       }
453f93e5bdSPeter Brune     }
469566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(pc->pmat, i, &ncols, &cols, &vals));
473f93e5bdSPeter Brune   }
4869b97631SPeter Brune   if (jac->symmetric) {
4969b97631SPeter Brune     for (i = xe - 1; i >= xs; i--) {
509566063dSJacob Faibussowitsch       PetscCall(MatGetRow(pc->pmat, i, &ncols, &cols, &vals));
5169b97631SPeter Brune       r    = xarray[i - xs];
5269b97631SPeter Brune       anrm = 0.;
5369b97631SPeter Brune       for (j = 0; j < ncols; j++) {
54ad540459SPierre Jolivet         if (cols[j] >= ys && cols[j] < ye) r -= yarray[cols[j] - ys] * vals[j];
5569b97631SPeter Brune         anrm += PetscRealPart(PetscSqr(vals[j]));
5669b97631SPeter Brune       }
5769b97631SPeter Brune       if (anrm > 0.) {
5869b97631SPeter Brune         for (j = 0; j < ncols; j++) {
59ad540459SPierre Jolivet           if (cols[j] >= ys && cols[j] < ye) yarray[cols[j] - ys] += vals[j] * lambda * r / anrm;
6069b97631SPeter Brune         }
6169b97631SPeter Brune       }
629566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(pc->pmat, i, &ncols, &cols, &vals));
6369b97631SPeter Brune     }
6469b97631SPeter Brune   }
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
683f93e5bdSPeter Brune }
693f93e5bdSPeter Brune 
PCSetFromOptions_Kaczmarz(PC pc,PetscOptionItems PetscOptionsObject)70ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_Kaczmarz(PC pc, PetscOptionItems PetscOptionsObject)
71d71ae5a4SJacob Faibussowitsch {
723f93e5bdSPeter Brune   PC_Kaczmarz *jac = (PC_Kaczmarz *)pc->data;
733f93e5bdSPeter Brune 
743f93e5bdSPeter Brune   PetscFunctionBegin;
75d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Kaczmarz options");
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_kaczmarz_lambda", "relaxation factor (0 < lambda)", "", jac->lambda, &jac->lambda, NULL));
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_kaczmarz_symmetric", "apply row projections symmetrically", "", jac->symmetric, &jac->symmetric, NULL));
78d0609cedSBarry Smith   PetscOptionsHeadEnd();
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
803f93e5bdSPeter Brune }
813f93e5bdSPeter Brune 
PCView_Kaczmarz(PC pc,PetscViewer viewer)8266976f2fSJacob Faibussowitsch static PetscErrorCode PCView_Kaczmarz(PC pc, PetscViewer viewer)
83d71ae5a4SJacob Faibussowitsch {
843f93e5bdSPeter Brune   PC_Kaczmarz *jac = (PC_Kaczmarz *)pc->data;
85*9f196a02SMartin Diehl   PetscBool    isascii;
863f93e5bdSPeter Brune 
873f93e5bdSPeter Brune   PetscFunctionBegin;
88*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
89*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  lambda = %g\n", (double)jac->lambda));
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
913f93e5bdSPeter Brune }
923f93e5bdSPeter Brune 
933f93e5bdSPeter Brune /*MC
941d27aa22SBarry Smith      PCKACZMARZ - Kaczmarz iteration {cite}`kaczmarz1937angenaherte`
953f93e5bdSPeter Brune 
96350f1385SJose E. Roman    Options Database Keys:
971d27aa22SBarry Smith +  -pc_kaczmarz_lambda <lambda> - Sets damping parameter defaults to 1.0
981d27aa22SBarry Smith -  -pc_kaczmarz_symmetric       - Apply the row projections symmetrically
993f93e5bdSPeter Brune 
1003f93e5bdSPeter Brune    Level: beginner
1013f93e5bdSPeter Brune 
102f1580f4eSBarry Smith    Note:
1031d27aa22SBarry Smith    In parallel this is block-Jacobi with Kaczmarz inner solve on each processor.
1043f93e5bdSPeter Brune 
105562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, `PCBJACOBI`
1063f93e5bdSPeter Brune M*/
1073f93e5bdSPeter Brune 
PCCreate_Kaczmarz(PC pc)108d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Kaczmarz(PC pc)
109d71ae5a4SJacob Faibussowitsch {
1103f93e5bdSPeter Brune   PC_Kaczmarz *jac;
1113f93e5bdSPeter Brune 
1123f93e5bdSPeter Brune   PetscFunctionBegin;
1134dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
1143f93e5bdSPeter Brune 
1153f93e5bdSPeter Brune   pc->ops->apply          = PCApply_Kaczmarz;
1163f93e5bdSPeter Brune   pc->ops->setfromoptions = PCSetFromOptions_Kaczmarz;
1170a545947SLisandro Dalcin   pc->ops->setup          = NULL;
1183f93e5bdSPeter Brune   pc->ops->view           = PCView_Kaczmarz;
1193f93e5bdSPeter Brune   pc->ops->destroy        = PCDestroy_Kaczmarz;
1203f93e5bdSPeter Brune   pc->data                = (void *)jac;
1213f93e5bdSPeter Brune   jac->lambda             = 1.0;
12269b97631SPeter Brune   jac->symmetric          = PETSC_FALSE;
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1243f93e5bdSPeter Brune }
125