xref: /petsc/src/snes/impls/gs/gssecant.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h>
2ef107cf5SPeter Brune 
SNESComputeNGSDefaultSecant(SNES snes,Vec X,Vec B,PetscCtx ctx)3*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGSDefaultSecant(SNES snes, Vec X, Vec B, PetscCtx ctx)
4d71ae5a4SJacob Faibussowitsch {
5be95d8f1SBarry Smith   SNES_NGS          *gs = (SNES_NGS *)snes->data;
6ef107cf5SPeter Brune   PetscInt           i, j, k, ncolors;
7ef107cf5SPeter Brune   DM                 dm;
8ef107cf5SPeter Brune   PetscBool          flg;
98a86d3c5SBarry Smith   ISColoring         coloring = gs->coloring;
10ef107cf5SPeter Brune   MatColoring        mc;
11ef107cf5SPeter Brune   Vec                W, G, F;
12737a7e12SPeter Brune   PetscScalar        h = gs->h;
13ef107cf5SPeter Brune   IS                *coloris;
14ef107cf5SPeter Brune   PetscScalar        f, g, x, w, d;
15c1132020SPeter Brune   PetscReal          dxt, xt, ft, ft1 = 0;
16ef107cf5SPeter Brune   const PetscInt    *idx;
17c03df580SPeter Brune   PetscInt           size, s;
18ef107cf5SPeter Brune   PetscReal          atol, rtol, stol;
19ef107cf5SPeter Brune   PetscInt           its;
206b72add0SBarry Smith   SNESFunctionFn    *func;
21ef107cf5SPeter Brune   void              *fctx;
22c03df580SPeter Brune   PetscBool          mat = gs->secant_mat, equal, isdone, alldone;
2378e24b28SBarry Smith   PetscScalar       *xa, *wa;
2478e24b28SBarry Smith   const PetscScalar *fa, *ga;
25ef107cf5SPeter Brune 
26ef107cf5SPeter Brune   PetscFunctionBegin;
2748a46eb9SPierre Jolivet   if (snes->nwork < 3) PetscCall(SNESSetWorkVecs(snes, 3));
28ef107cf5SPeter Brune   W = snes->work[0];
29ef107cf5SPeter Brune   G = snes->work[1];
30ef107cf5SPeter Brune   F = snes->work[2];
319566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &s, NULL));
329566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
339566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
349566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, NULL, &func, &fctx));
358a86d3c5SBarry Smith   if (!coloring) {
36ef107cf5SPeter Brune     /* create the coloring */
379566063dSJacob Faibussowitsch     PetscCall(DMHasColoring(dm, &flg));
38737a7e12SPeter Brune     if (flg && !mat) {
399566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, IS_COLORING_GLOBAL, &coloring));
40ef107cf5SPeter Brune     } else {
419566063dSJacob Faibussowitsch       if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));
429566063dSJacob Faibussowitsch       PetscCall(MatColoringCreate(snes->jacobian, &mc));
439566063dSJacob Faibussowitsch       PetscCall(MatColoringSetDistance(mc, 1));
449566063dSJacob Faibussowitsch       PetscCall(MatColoringSetFromOptions(mc));
459566063dSJacob Faibussowitsch       PetscCall(MatColoringApply(mc, &coloring));
469566063dSJacob Faibussowitsch       PetscCall(MatColoringDestroy(&mc));
47ef107cf5SPeter Brune     }
488a86d3c5SBarry Smith     gs->coloring = coloring;
49ef107cf5SPeter Brune   }
509566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(coloring, PETSC_USE_POINTER, &ncolors, &coloris));
519566063dSJacob Faibussowitsch   PetscCall(VecEqual(X, snes->vec_sol, &equal));
52c03df580SPeter Brune   if (equal && snes->normschedule == SNES_NORM_ALWAYS) {
53737a7e12SPeter Brune     /* assume that the function is already computed */
549566063dSJacob Faibussowitsch     PetscCall(VecCopy(snes->vec_func, F));
55737a7e12SPeter Brune   } else {
569566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_NGSFuncEval, snes, X, B, 0));
579566063dSJacob Faibussowitsch     PetscCall((*func)(snes, X, F, fctx));
589566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_NGSFuncEval, snes, X, B, 0));
599566063dSJacob Faibussowitsch     if (B) PetscCall(VecAXPY(F, -1.0, B));
60737a7e12SPeter Brune   }
61737a7e12SPeter Brune   for (i = 0; i < ncolors; i++) {
629566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(coloris[i], &idx));
639566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(coloris[i], &size));
649566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, W));
659566063dSJacob Faibussowitsch     PetscCall(VecGetArray(W, &wa));
66ad540459SPierre Jolivet     for (j = 0; j < size; j++) wa[idx[j] - s] += h;
679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(W, &wa));
689566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_NGSFuncEval, snes, X, B, 0));
699566063dSJacob Faibussowitsch     PetscCall((*func)(snes, W, G, fctx));
709566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_NGSFuncEval, snes, X, B, 0));
719566063dSJacob Faibussowitsch     if (B) PetscCall(VecAXPY(G, -1.0, B));
72737a7e12SPeter Brune     for (k = 0; k < its; k++) {
73737a7e12SPeter Brune       dxt = 0.;
74737a7e12SPeter Brune       xt  = 0.;
75c03df580SPeter Brune       ft  = 0.;
769566063dSJacob Faibussowitsch       PetscCall(VecGetArray(W, &wa));
779566063dSJacob Faibussowitsch       PetscCall(VecGetArray(X, &xa));
789566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(F, &fa));
799566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(G, &ga));
80ef107cf5SPeter Brune       for (j = 0; j < size; j++) {
81c03df580SPeter Brune         f = fa[idx[j] - s];
82c03df580SPeter Brune         x = xa[idx[j] - s];
83c03df580SPeter Brune         g = ga[idx[j] - s];
84c03df580SPeter Brune         w = wa[idx[j] - s];
85737a7e12SPeter Brune         if (PetscAbsScalar(g - f) > atol) {
864b655d89SMatthew G. Knepley           /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */
87ef107cf5SPeter Brune           d = (x * g - w * f) / PetscRealPart(g - f);
88737a7e12SPeter Brune         } else {
89737a7e12SPeter Brune           d = x;
90737a7e12SPeter Brune         }
91400ea6b8SPeter Brune         dxt += PetscRealPart(PetscSqr(d - x));
92400ea6b8SPeter Brune         xt += PetscRealPart(PetscSqr(x));
93400ea6b8SPeter Brune         ft += PetscRealPart(PetscSqr(f));
94c03df580SPeter Brune         xa[idx[j] - s] = d;
95ef107cf5SPeter Brune       }
969566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(X, &xa));
979566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(F, &fa));
989566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(G, &ga));
999566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(W, &wa));
100c03df580SPeter Brune 
101315bab25SPeter Brune       if (k == 0) ft1 = PetscSqrtReal(ft);
102c03df580SPeter Brune       if (k < its - 1) {
103c03df580SPeter Brune         isdone = PETSC_FALSE;
104c03df580SPeter Brune         if (stol * PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE;
105c03df580SPeter Brune         if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE;
106c03df580SPeter Brune         if (rtol * ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE;
1075440e5dcSBarry Smith         PetscCallMPI(MPIU_Allreduce(&isdone, &alldone, 1, MPI_C_BOOL, MPI_LAND, PetscObjectComm((PetscObject)snes)));
108c03df580SPeter Brune         if (alldone) break;
109c03df580SPeter Brune       }
110737a7e12SPeter Brune       if (i < ncolors - 1 || k < its - 1) {
1119566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NGSFuncEval, snes, X, B, 0));
1129566063dSJacob Faibussowitsch         PetscCall((*func)(snes, X, F, fctx));
1139566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NGSFuncEval, snes, X, B, 0));
1149566063dSJacob Faibussowitsch         if (B) PetscCall(VecAXPY(F, -1.0, B));
115737a7e12SPeter Brune       }
116737a7e12SPeter Brune       if (k < its - 1) {
1179566063dSJacob Faibussowitsch         PetscCall(VecSwap(X, W));
1189566063dSJacob Faibussowitsch         PetscCall(VecSwap(F, G));
119737a7e12SPeter Brune       }
120ef107cf5SPeter Brune     }
121ef107cf5SPeter Brune   }
1229566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(coloring, PETSC_USE_POINTER, &coloris));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124ef107cf5SPeter Brune }
125