1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> 2ef107cf5SPeter Brune 3ef107cf5SPeter Brune #undef __FUNCT__ 4235fd6e6SBarry Smith #define __FUNCT__ "SNESNGSDestroy_Private" 5235fd6e6SBarry Smith static PetscErrorCode SNESNGSDestroy_Private(ISColoring coloring) 6ef107cf5SPeter Brune { 7ef107cf5SPeter Brune PetscErrorCode ierr; 8ef107cf5SPeter Brune 9ef107cf5SPeter Brune PetscFunctionBegin; 10ef107cf5SPeter Brune ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 11ef107cf5SPeter Brune PetscFunctionReturn(0); 12ef107cf5SPeter Brune } 13ef107cf5SPeter Brune 14ef107cf5SPeter Brune #undef __FUNCT__ 15be95d8f1SBarry Smith #define __FUNCT__ "SNESComputeNGSDefaultSecant" 16be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx) 17ef107cf5SPeter Brune { 18ef107cf5SPeter Brune PetscErrorCode ierr; 19be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 20ef107cf5SPeter Brune PetscInt i,j,k,ncolors; 21ef107cf5SPeter Brune DM dm; 22ef107cf5SPeter Brune PetscBool flg; 23ef107cf5SPeter Brune ISColoring coloring; 24ef107cf5SPeter Brune MatColoring mc; 25ef107cf5SPeter Brune Vec W,G,F; 26737a7e12SPeter Brune PetscScalar h=gs->h; 27ef107cf5SPeter Brune IS *coloris; 28ef107cf5SPeter Brune PetscScalar f,g,x,w,d; 29c1132020SPeter Brune PetscReal dxt,xt,ft,ft1=0; 30ef107cf5SPeter Brune const PetscInt *idx; 31c03df580SPeter Brune PetscInt size,s; 32ef107cf5SPeter Brune PetscReal atol,rtol,stol; 33ef107cf5SPeter Brune PetscInt its; 34ef107cf5SPeter Brune PetscErrorCode (*func)(SNES,Vec,Vec,void*); 35ef107cf5SPeter Brune void *fctx; 36ef107cf5SPeter Brune PetscContainer colorcontainer; 37c03df580SPeter Brune PetscBool mat = gs->secant_mat,equal,isdone,alldone; 38c03df580SPeter Brune PetscScalar *xa,*fa,*wa,*ga; 39ef107cf5SPeter Brune 40ef107cf5SPeter Brune PetscFunctionBegin; 41ef107cf5SPeter Brune if (snes->nwork < 3) { 42ef107cf5SPeter Brune ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 43ef107cf5SPeter Brune } 44ef107cf5SPeter Brune W = snes->work[0]; 45ef107cf5SPeter Brune G = snes->work[1]; 46ef107cf5SPeter Brune F = snes->work[2]; 47c03df580SPeter Brune ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr); 48be95d8f1SBarry Smith ierr = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 49ef107cf5SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 50ef107cf5SPeter Brune ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr); 51be95d8f1SBarry Smith ierr = PetscObjectQuery((PetscObject)snes,"SNESNGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr); 52ef107cf5SPeter Brune if (!colorcontainer) { 53ef107cf5SPeter Brune /* create the coloring */ 54ef107cf5SPeter Brune ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr); 55737a7e12SPeter Brune if (flg && !mat) { 56ef107cf5SPeter Brune ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr); 57ef107cf5SPeter Brune } else { 58ef107cf5SPeter Brune if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);} 59ef107cf5SPeter Brune ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr); 60ef107cf5SPeter Brune ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr); 61ef107cf5SPeter Brune ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 62ef107cf5SPeter Brune ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr); 63ef107cf5SPeter Brune ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 64ef107cf5SPeter Brune } 65ef107cf5SPeter Brune ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr); 66ef107cf5SPeter Brune ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr); 67235fd6e6SBarry Smith ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))SNESNGSDestroy_Private);CHKERRQ(ierr); 68be95d8f1SBarry Smith ierr = PetscObjectCompose((PetscObject)snes,"SNESNGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr); 69ef107cf5SPeter Brune ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr); 70ef107cf5SPeter Brune } else { 71ef107cf5SPeter Brune ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr); 72ef107cf5SPeter Brune } 73ef107cf5SPeter Brune ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr); 74737a7e12SPeter Brune ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr); 75c03df580SPeter Brune if (equal && snes->normschedule == SNES_NORM_ALWAYS) { 76737a7e12SPeter Brune /* assume that the function is already computed */ 77737a7e12SPeter Brune ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr); 78737a7e12SPeter Brune } else { 79be95d8f1SBarry Smith ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 80ef107cf5SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 81be95d8f1SBarry Smith ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 82ef107cf5SPeter Brune if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 83737a7e12SPeter Brune } 84c03df580SPeter Brune ierr = VecGetArray(X,&xa);CHKERRQ(ierr); 85c03df580SPeter Brune ierr = VecGetArray(F,&fa);CHKERRQ(ierr); 86c03df580SPeter Brune ierr = VecGetArray(G,&ga);CHKERRQ(ierr); 87c03df580SPeter Brune ierr = VecGetArray(W,&wa);CHKERRQ(ierr); 88737a7e12SPeter Brune for (i=0;i<ncolors;i++) { 89ef107cf5SPeter Brune ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr); 90ef107cf5SPeter Brune ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr); 91ef107cf5SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 92ef107cf5SPeter Brune for (j=0;j<size;j++) { 93c03df580SPeter Brune wa[idx[j]-s] += h; 94ef107cf5SPeter Brune } 95be95d8f1SBarry Smith ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 96ef107cf5SPeter Brune ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr); 97be95d8f1SBarry Smith ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 98ef107cf5SPeter Brune if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);} 99737a7e12SPeter Brune for (k=0;k<its;k++) { 100737a7e12SPeter Brune dxt = 0.; 101737a7e12SPeter Brune xt = 0.; 102c03df580SPeter Brune ft = 0.; 103ef107cf5SPeter Brune for (j=0;j<size;j++) { 104c03df580SPeter Brune f = fa[idx[j]-s]; 105c03df580SPeter Brune x = xa[idx[j]-s]; 106c03df580SPeter Brune g = ga[idx[j]-s]; 107c03df580SPeter Brune w = wa[idx[j]-s]; 108737a7e12SPeter Brune if (PetscAbsScalar(g-f) > atol) { 109*4b655d89SMatthew G. Knepley /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */ 110ef107cf5SPeter Brune d = (x*g-w*f) / PetscRealPart(g-f); 111737a7e12SPeter Brune } else { 112737a7e12SPeter Brune d = x; 113737a7e12SPeter Brune } 114400ea6b8SPeter Brune dxt += PetscRealPart(PetscSqr(d-x)); 115400ea6b8SPeter Brune xt += PetscRealPart(PetscSqr(x)); 116400ea6b8SPeter Brune ft += PetscRealPart(PetscSqr(f)); 117c03df580SPeter Brune xa[idx[j]-s] = d; 118ef107cf5SPeter Brune } 119c03df580SPeter Brune 120315bab25SPeter Brune if (k == 0) ft1 = PetscSqrtReal(ft); 121c03df580SPeter Brune if (k<its-1) { 122c03df580SPeter Brune isdone = PETSC_FALSE; 123c03df580SPeter Brune if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE; 124c03df580SPeter Brune if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE; 125c03df580SPeter Brune if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE; 126c03df580SPeter Brune ierr = MPI_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 127c03df580SPeter Brune if (alldone) break; 128c03df580SPeter Brune } 129737a7e12SPeter Brune if (i < ncolors-1 || k < its-1) { 130be95d8f1SBarry Smith ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 131737a7e12SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 132be95d8f1SBarry Smith ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr); 133737a7e12SPeter Brune if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 134737a7e12SPeter Brune } 135737a7e12SPeter Brune if (k<its-1) { 136737a7e12SPeter Brune ierr = VecSwap(X,W);CHKERRQ(ierr); 137737a7e12SPeter Brune ierr = VecSwap(F,G);CHKERRQ(ierr); 138737a7e12SPeter Brune } 139ef107cf5SPeter Brune } 140ef107cf5SPeter Brune } 141c03df580SPeter Brune ierr = VecRestoreArray(X,&xa);CHKERRQ(ierr); 142c03df580SPeter Brune ierr = VecRestoreArray(F,&fa);CHKERRQ(ierr); 143c03df580SPeter Brune ierr = VecRestoreArray(G,&ga);CHKERRQ(ierr); 144c03df580SPeter Brune ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr); 145ef107cf5SPeter Brune ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr); 146ef107cf5SPeter Brune PetscFunctionReturn(0); 147ef107cf5SPeter Brune } 148