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 { 5ef107cf5SPeter Brune PetscErrorCode ierr; 6be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 7ef107cf5SPeter Brune PetscInt i,j,k,ncolors; 8ef107cf5SPeter Brune DM dm; 9ef107cf5SPeter Brune PetscBool flg; 108a86d3c5SBarry Smith ISColoring coloring = gs->coloring; 11ef107cf5SPeter Brune MatColoring mc; 12ef107cf5SPeter Brune Vec W,G,F; 13737a7e12SPeter Brune PetscScalar h=gs->h; 14ef107cf5SPeter Brune IS *coloris; 15ef107cf5SPeter Brune PetscScalar f,g,x,w,d; 16c1132020SPeter Brune PetscReal dxt,xt,ft,ft1=0; 17ef107cf5SPeter Brune const PetscInt *idx; 18c03df580SPeter Brune PetscInt size,s; 19ef107cf5SPeter Brune PetscReal atol,rtol,stol; 20ef107cf5SPeter Brune PetscInt its; 21ef107cf5SPeter Brune PetscErrorCode (*func)(SNES,Vec,Vec,void*); 22ef107cf5SPeter Brune void *fctx; 23c03df580SPeter Brune PetscBool mat = gs->secant_mat,equal,isdone,alldone; 24*78e24b28SBarry Smith PetscScalar *xa,*wa; 25*78e24b28SBarry Smith const PetscScalar *fa,*ga; 26ef107cf5SPeter Brune 27ef107cf5SPeter Brune PetscFunctionBegin; 28ef107cf5SPeter Brune if (snes->nwork < 3) { 29ef107cf5SPeter Brune ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 30ef107cf5SPeter Brune } 31ef107cf5SPeter Brune W = snes->work[0]; 32ef107cf5SPeter Brune G = snes->work[1]; 33ef107cf5SPeter Brune F = snes->work[2]; 34c03df580SPeter Brune ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr); 35be95d8f1SBarry Smith ierr = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 36ef107cf5SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 37ef107cf5SPeter Brune ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr); 388a86d3c5SBarry Smith if (!coloring) { 39ef107cf5SPeter Brune /* create the coloring */ 40ef107cf5SPeter Brune ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr); 41737a7e12SPeter Brune if (flg && !mat) { 42ef107cf5SPeter Brune ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr); 43ef107cf5SPeter Brune } else { 44ef107cf5SPeter Brune if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);} 45ef107cf5SPeter Brune ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr); 46ef107cf5SPeter Brune ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr); 47ef107cf5SPeter Brune ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 48ef107cf5SPeter Brune ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr); 49ef107cf5SPeter Brune ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 50ef107cf5SPeter Brune } 518a86d3c5SBarry Smith gs->coloring = coloring; 52ef107cf5SPeter Brune } 53071fcb05SBarry Smith ierr = ISColoringGetIS(coloring,PETSC_USE_POINTER,&ncolors,&coloris);CHKERRQ(ierr); 54737a7e12SPeter Brune ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr); 55c03df580SPeter Brune if (equal && snes->normschedule == SNES_NORM_ALWAYS) { 56737a7e12SPeter Brune /* assume that the function is already computed */ 57737a7e12SPeter Brune ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr); 58737a7e12SPeter Brune } else { 59be95d8f1SBarry Smith ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 60ef107cf5SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 61be95d8f1SBarry Smith ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 62ef107cf5SPeter Brune if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 63737a7e12SPeter Brune } 64737a7e12SPeter Brune for (i=0;i<ncolors;i++) { 65ef107cf5SPeter Brune ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr); 66ef107cf5SPeter Brune ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr); 67ef107cf5SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 68*78e24b28SBarry Smith ierr = VecGetArray(W,&wa);CHKERRQ(ierr); 69ef107cf5SPeter Brune for (j=0;j<size;j++) { 70c03df580SPeter Brune wa[idx[j]-s] += h; 71ef107cf5SPeter Brune } 72*78e24b28SBarry Smith ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr); 73be95d8f1SBarry Smith ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 74ef107cf5SPeter Brune ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr); 75be95d8f1SBarry Smith ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 76ef107cf5SPeter Brune if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);} 77737a7e12SPeter Brune for (k=0;k<its;k++) { 78737a7e12SPeter Brune dxt = 0.; 79737a7e12SPeter Brune xt = 0.; 80c03df580SPeter Brune ft = 0.; 81*78e24b28SBarry Smith ierr = VecGetArray(W,&wa);CHKERRQ(ierr); 82*78e24b28SBarry Smith ierr = VecGetArray(X,&xa);CHKERRQ(ierr); 83*78e24b28SBarry Smith ierr = VecGetArrayRead(F,&fa);CHKERRQ(ierr); 84*78e24b28SBarry Smith ierr = VecGetArrayRead(G,&ga);CHKERRQ(ierr); 85ef107cf5SPeter Brune for (j=0;j<size;j++) { 86c03df580SPeter Brune f = fa[idx[j]-s]; 87c03df580SPeter Brune x = xa[idx[j]-s]; 88c03df580SPeter Brune g = ga[idx[j]-s]; 89c03df580SPeter Brune w = wa[idx[j]-s]; 90737a7e12SPeter Brune if (PetscAbsScalar(g-f) > atol) { 914b655d89SMatthew G. Knepley /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */ 92ef107cf5SPeter Brune d = (x*g-w*f) / PetscRealPart(g-f); 93737a7e12SPeter Brune } else { 94737a7e12SPeter Brune d = x; 95737a7e12SPeter Brune } 96400ea6b8SPeter Brune dxt += PetscRealPart(PetscSqr(d-x)); 97400ea6b8SPeter Brune xt += PetscRealPart(PetscSqr(x)); 98400ea6b8SPeter Brune ft += PetscRealPart(PetscSqr(f)); 99c03df580SPeter Brune xa[idx[j]-s] = d; 100ef107cf5SPeter Brune } 101*78e24b28SBarry Smith ierr = VecRestoreArray(X,&xa);CHKERRQ(ierr); 102*78e24b28SBarry Smith ierr = VecRestoreArrayRead(F,&fa);CHKERRQ(ierr); 103*78e24b28SBarry Smith ierr = VecRestoreArrayRead(G,&ga);CHKERRQ(ierr); 104*78e24b28SBarry Smith ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr); 105c03df580SPeter Brune 106315bab25SPeter Brune if (k == 0) ft1 = PetscSqrtReal(ft); 107c03df580SPeter Brune if (k<its-1) { 108c03df580SPeter Brune isdone = PETSC_FALSE; 109c03df580SPeter Brune if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE; 110c03df580SPeter Brune if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE; 111c03df580SPeter Brune if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE; 112b2566f29SBarry Smith ierr = MPIU_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 113c03df580SPeter Brune if (alldone) break; 114c03df580SPeter Brune } 115737a7e12SPeter Brune if (i < ncolors-1 || k < its-1) { 116be95d8f1SBarry Smith ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 117737a7e12SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 118be95d8f1SBarry Smith ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 119737a7e12SPeter Brune if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 120737a7e12SPeter Brune } 121737a7e12SPeter Brune if (k<its-1) { 122737a7e12SPeter Brune ierr = VecSwap(X,W);CHKERRQ(ierr); 123737a7e12SPeter Brune ierr = VecSwap(F,G);CHKERRQ(ierr); 124737a7e12SPeter Brune } 125ef107cf5SPeter Brune } 126ef107cf5SPeter Brune } 127071fcb05SBarry Smith ierr = ISColoringRestoreIS(coloring,PETSC_USE_POINTER,&coloris);CHKERRQ(ierr); 128ef107cf5SPeter Brune PetscFunctionReturn(0); 129ef107cf5SPeter Brune } 130