1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> 2ef107cf5SPeter Brune 39371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESComputeNGSDefaultSecant(SNES snes, Vec X, Vec B, void *ctx) { 4be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 5ef107cf5SPeter Brune PetscInt i, j, k, ncolors; 6ef107cf5SPeter Brune DM dm; 7ef107cf5SPeter Brune PetscBool flg; 88a86d3c5SBarry Smith ISColoring coloring = gs->coloring; 9ef107cf5SPeter Brune MatColoring mc; 10ef107cf5SPeter Brune Vec W, G, F; 11737a7e12SPeter Brune PetscScalar h = gs->h; 12ef107cf5SPeter Brune IS *coloris; 13ef107cf5SPeter Brune PetscScalar f, g, x, w, d; 14c1132020SPeter Brune PetscReal dxt, xt, ft, ft1 = 0; 15ef107cf5SPeter Brune const PetscInt *idx; 16c03df580SPeter Brune PetscInt size, s; 17ef107cf5SPeter Brune PetscReal atol, rtol, stol; 18ef107cf5SPeter Brune PetscInt its; 19ef107cf5SPeter Brune PetscErrorCode (*func)(SNES, Vec, Vec, void *); 20ef107cf5SPeter Brune void *fctx; 21c03df580SPeter Brune PetscBool mat = gs->secant_mat, equal, isdone, alldone; 2278e24b28SBarry Smith PetscScalar *xa, *wa; 2378e24b28SBarry Smith const PetscScalar *fa, *ga; 24ef107cf5SPeter Brune 25ef107cf5SPeter Brune PetscFunctionBegin; 26*48a46eb9SPierre Jolivet if (snes->nwork < 3) PetscCall(SNESSetWorkVecs(snes, 3)); 27ef107cf5SPeter Brune W = snes->work[0]; 28ef107cf5SPeter Brune G = snes->work[1]; 29ef107cf5SPeter Brune F = snes->work[2]; 309566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &s, NULL)); 319566063dSJacob Faibussowitsch PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its)); 329566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 339566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, NULL, &func, &fctx)); 348a86d3c5SBarry Smith if (!coloring) { 35ef107cf5SPeter Brune /* create the coloring */ 369566063dSJacob Faibussowitsch PetscCall(DMHasColoring(dm, &flg)); 37737a7e12SPeter Brune if (flg && !mat) { 389566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm, IS_COLORING_GLOBAL, &coloring)); 39ef107cf5SPeter Brune } else { 409566063dSJacob Faibussowitsch if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes)); 419566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(snes->jacobian, &mc)); 429566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(mc, 1)); 439566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 449566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc, &coloring)); 459566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 46ef107cf5SPeter Brune } 478a86d3c5SBarry Smith gs->coloring = coloring; 48ef107cf5SPeter Brune } 499566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(coloring, PETSC_USE_POINTER, &ncolors, &coloris)); 509566063dSJacob Faibussowitsch PetscCall(VecEqual(X, snes->vec_sol, &equal)); 51c03df580SPeter Brune if (equal && snes->normschedule == SNES_NORM_ALWAYS) { 52737a7e12SPeter Brune /* assume that the function is already computed */ 539566063dSJacob Faibussowitsch PetscCall(VecCopy(snes->vec_func, F)); 54737a7e12SPeter Brune } else { 559566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NGSFuncEval, snes, X, B, 0)); 569566063dSJacob Faibussowitsch PetscCall((*func)(snes, X, F, fctx)); 579566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NGSFuncEval, snes, X, B, 0)); 589566063dSJacob Faibussowitsch if (B) PetscCall(VecAXPY(F, -1.0, B)); 59737a7e12SPeter Brune } 60737a7e12SPeter Brune for (i = 0; i < ncolors; i++) { 619566063dSJacob Faibussowitsch PetscCall(ISGetIndices(coloris[i], &idx)); 629566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(coloris[i], &size)); 639566063dSJacob Faibussowitsch PetscCall(VecCopy(X, W)); 649566063dSJacob Faibussowitsch PetscCall(VecGetArray(W, &wa)); 659371c9d4SSatish Balay for (j = 0; j < size; j++) { wa[idx[j] - s] += h; } 669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(W, &wa)); 679566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NGSFuncEval, snes, X, B, 0)); 689566063dSJacob Faibussowitsch PetscCall((*func)(snes, W, G, fctx)); 699566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NGSFuncEval, snes, X, B, 0)); 709566063dSJacob Faibussowitsch if (B) PetscCall(VecAXPY(G, -1.0, B)); 71737a7e12SPeter Brune for (k = 0; k < its; k++) { 72737a7e12SPeter Brune dxt = 0.; 73737a7e12SPeter Brune xt = 0.; 74c03df580SPeter Brune ft = 0.; 759566063dSJacob Faibussowitsch PetscCall(VecGetArray(W, &wa)); 769566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &xa)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &fa)); 789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(G, &ga)); 79ef107cf5SPeter Brune for (j = 0; j < size; j++) { 80c03df580SPeter Brune f = fa[idx[j] - s]; 81c03df580SPeter Brune x = xa[idx[j] - s]; 82c03df580SPeter Brune g = ga[idx[j] - s]; 83c03df580SPeter Brune w = wa[idx[j] - s]; 84737a7e12SPeter Brune if (PetscAbsScalar(g - f) > atol) { 854b655d89SMatthew G. Knepley /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */ 86ef107cf5SPeter Brune d = (x * g - w * f) / PetscRealPart(g - f); 87737a7e12SPeter Brune } else { 88737a7e12SPeter Brune d = x; 89737a7e12SPeter Brune } 90400ea6b8SPeter Brune dxt += PetscRealPart(PetscSqr(d - x)); 91400ea6b8SPeter Brune xt += PetscRealPart(PetscSqr(x)); 92400ea6b8SPeter Brune ft += PetscRealPart(PetscSqr(f)); 93c03df580SPeter Brune xa[idx[j] - s] = d; 94ef107cf5SPeter Brune } 959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &xa)); 969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &fa)); 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(G, &ga)); 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(W, &wa)); 99c03df580SPeter Brune 100315bab25SPeter Brune if (k == 0) ft1 = PetscSqrtReal(ft); 101c03df580SPeter Brune if (k < its - 1) { 102c03df580SPeter Brune isdone = PETSC_FALSE; 103c03df580SPeter Brune if (stol * PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE; 104c03df580SPeter Brune if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE; 105c03df580SPeter Brune if (rtol * ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE; 1061c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&isdone, &alldone, 1, MPIU_BOOL, MPI_BAND, PetscObjectComm((PetscObject)snes))); 107c03df580SPeter Brune if (alldone) break; 108c03df580SPeter Brune } 109737a7e12SPeter Brune if (i < ncolors - 1 || k < its - 1) { 1109566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NGSFuncEval, snes, X, B, 0)); 1119566063dSJacob Faibussowitsch PetscCall((*func)(snes, X, F, fctx)); 1129566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NGSFuncEval, snes, X, B, 0)); 1139566063dSJacob Faibussowitsch if (B) PetscCall(VecAXPY(F, -1.0, B)); 114737a7e12SPeter Brune } 115737a7e12SPeter Brune if (k < its - 1) { 1169566063dSJacob Faibussowitsch PetscCall(VecSwap(X, W)); 1179566063dSJacob Faibussowitsch PetscCall(VecSwap(F, G)); 118737a7e12SPeter Brune } 119ef107cf5SPeter Brune } 120ef107cf5SPeter Brune } 1219566063dSJacob Faibussowitsch PetscCall(ISColoringRestoreIS(coloring, PETSC_USE_POINTER, &coloris)); 122ef107cf5SPeter Brune PetscFunctionReturn(0); 123ef107cf5SPeter Brune } 124