1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> 2ef107cf5SPeter Brune 3be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx) 4ef107cf5SPeter Brune { 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; 20ef107cf5SPeter Brune PetscErrorCode (*func)(SNES,Vec,Vec,void*); 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; 27ef107cf5SPeter Brune if (snes->nwork < 3) { 28*9566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes,3)); 29ef107cf5SPeter Brune } 30ef107cf5SPeter Brune W = snes->work[0]; 31ef107cf5SPeter Brune G = snes->work[1]; 32ef107cf5SPeter Brune F = snes->work[2]; 33*9566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X,&s,NULL)); 34*9566063dSJacob Faibussowitsch PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its)); 35*9566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 36*9566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,NULL,&func,&fctx)); 378a86d3c5SBarry Smith if (!coloring) { 38ef107cf5SPeter Brune /* create the coloring */ 39*9566063dSJacob Faibussowitsch PetscCall(DMHasColoring(dm,&flg)); 40737a7e12SPeter Brune if (flg && !mat) { 41*9566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring)); 42ef107cf5SPeter Brune } else { 43*9566063dSJacob Faibussowitsch if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes)); 44*9566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(snes->jacobian,&mc)); 45*9566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(mc,1)); 46*9566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 47*9566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc,&coloring)); 48*9566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 49ef107cf5SPeter Brune } 508a86d3c5SBarry Smith gs->coloring = coloring; 51ef107cf5SPeter Brune } 52*9566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(coloring,PETSC_USE_POINTER,&ncolors,&coloris)); 53*9566063dSJacob Faibussowitsch PetscCall(VecEqual(X,snes->vec_sol,&equal)); 54c03df580SPeter Brune if (equal && snes->normschedule == SNES_NORM_ALWAYS) { 55737a7e12SPeter Brune /* assume that the function is already computed */ 56*9566063dSJacob Faibussowitsch PetscCall(VecCopy(snes->vec_func,F)); 57737a7e12SPeter Brune } else { 58*9566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0)); 59*9566063dSJacob Faibussowitsch PetscCall((*func)(snes,X,F,fctx)); 60*9566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0)); 61*9566063dSJacob Faibussowitsch if (B) PetscCall(VecAXPY(F,-1.0,B)); 62737a7e12SPeter Brune } 63737a7e12SPeter Brune for (i=0;i<ncolors;i++) { 64*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(coloris[i],&idx)); 65*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(coloris[i],&size)); 66*9566063dSJacob Faibussowitsch PetscCall(VecCopy(X,W)); 67*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(W,&wa)); 68ef107cf5SPeter Brune for (j=0;j<size;j++) { 69c03df580SPeter Brune wa[idx[j]-s] += h; 70ef107cf5SPeter Brune } 71*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(W,&wa)); 72*9566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0)); 73*9566063dSJacob Faibussowitsch PetscCall((*func)(snes,W,G,fctx)); 74*9566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0)); 75*9566063dSJacob Faibussowitsch if (B) PetscCall(VecAXPY(G,-1.0,B)); 76737a7e12SPeter Brune for (k=0;k<its;k++) { 77737a7e12SPeter Brune dxt = 0.; 78737a7e12SPeter Brune xt = 0.; 79c03df580SPeter Brune ft = 0.; 80*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(W,&wa)); 81*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&xa)); 82*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F,&fa)); 83*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(G,&ga)); 84ef107cf5SPeter Brune for (j=0;j<size;j++) { 85c03df580SPeter Brune f = fa[idx[j]-s]; 86c03df580SPeter Brune x = xa[idx[j]-s]; 87c03df580SPeter Brune g = ga[idx[j]-s]; 88c03df580SPeter Brune w = wa[idx[j]-s]; 89737a7e12SPeter Brune if (PetscAbsScalar(g-f) > atol) { 904b655d89SMatthew G. Knepley /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */ 91ef107cf5SPeter Brune d = (x*g-w*f) / PetscRealPart(g-f); 92737a7e12SPeter Brune } else { 93737a7e12SPeter Brune d = x; 94737a7e12SPeter Brune } 95400ea6b8SPeter Brune dxt += PetscRealPart(PetscSqr(d-x)); 96400ea6b8SPeter Brune xt += PetscRealPart(PetscSqr(x)); 97400ea6b8SPeter Brune ft += PetscRealPart(PetscSqr(f)); 98c03df580SPeter Brune xa[idx[j]-s] = d; 99ef107cf5SPeter Brune } 100*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&xa)); 101*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F,&fa)); 102*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(G,&ga)); 103*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(W,&wa)); 104c03df580SPeter Brune 105315bab25SPeter Brune if (k == 0) ft1 = PetscSqrtReal(ft); 106c03df580SPeter Brune if (k<its-1) { 107c03df580SPeter Brune isdone = PETSC_FALSE; 108c03df580SPeter Brune if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE; 109c03df580SPeter Brune if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE; 110c03df580SPeter Brune if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE; 111*9566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes))); 112c03df580SPeter Brune if (alldone) break; 113c03df580SPeter Brune } 114737a7e12SPeter Brune if (i < ncolors-1 || k < its-1) { 115*9566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0)); 116*9566063dSJacob Faibussowitsch PetscCall((*func)(snes,X,F,fctx)); 117*9566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0)); 118*9566063dSJacob Faibussowitsch if (B) PetscCall(VecAXPY(F,-1.0,B)); 119737a7e12SPeter Brune } 120737a7e12SPeter Brune if (k<its-1) { 121*9566063dSJacob Faibussowitsch PetscCall(VecSwap(X,W)); 122*9566063dSJacob Faibussowitsch PetscCall(VecSwap(F,G)); 123737a7e12SPeter Brune } 124ef107cf5SPeter Brune } 125ef107cf5SPeter Brune } 126*9566063dSJacob Faibussowitsch PetscCall(ISColoringRestoreIS(coloring,PETSC_USE_POINTER,&coloris)); 127ef107cf5SPeter Brune PetscFunctionReturn(0); 128ef107cf5SPeter Brune } 129