1*ef107cf5SPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2*ef107cf5SPeter Brune #include <petscdm.h> 3*ef107cf5SPeter Brune 4*ef107cf5SPeter Brune #undef __FUNCT__ 5*ef107cf5SPeter Brune #define __FUNCT__ "GSDestroy_Private" 6*ef107cf5SPeter Brune PetscErrorCode GSDestroy_Private(ISColoring coloring) 7*ef107cf5SPeter Brune { 8*ef107cf5SPeter Brune PetscErrorCode ierr; 9*ef107cf5SPeter Brune 10*ef107cf5SPeter Brune PetscFunctionBegin; 11*ef107cf5SPeter Brune ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 12*ef107cf5SPeter Brune PetscFunctionReturn(0); 13*ef107cf5SPeter Brune 14*ef107cf5SPeter Brune } 15*ef107cf5SPeter Brune 16*ef107cf5SPeter Brune #undef __FUNCT__ 17*ef107cf5SPeter Brune #define __FUNCT__ "SNESComputeGSDefaultSecant" 18*ef107cf5SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx) 19*ef107cf5SPeter Brune { 20*ef107cf5SPeter Brune 21*ef107cf5SPeter Brune PetscErrorCode ierr; 22*ef107cf5SPeter Brune PetscInt i,j,k,ncolors; 23*ef107cf5SPeter Brune DM dm; 24*ef107cf5SPeter Brune PetscBool flg; 25*ef107cf5SPeter Brune ISColoring coloring; 26*ef107cf5SPeter Brune MatColoring mc; 27*ef107cf5SPeter Brune Vec W,G,F; 28*ef107cf5SPeter Brune PetscScalar h=1e-6; 29*ef107cf5SPeter Brune IS *coloris; 30*ef107cf5SPeter Brune PetscScalar f,g,x,w,d; 31*ef107cf5SPeter Brune const PetscInt *idx; 32*ef107cf5SPeter Brune PetscInt size; 33*ef107cf5SPeter Brune PetscReal atol,rtol,stol; 34*ef107cf5SPeter Brune PetscInt its; 35*ef107cf5SPeter Brune PetscErrorCode (*func)(SNES,Vec,Vec,void*); 36*ef107cf5SPeter Brune void *fctx; 37*ef107cf5SPeter Brune PetscContainer colorcontainer; 38*ef107cf5SPeter Brune 39*ef107cf5SPeter Brune PetscFunctionBegin; 40*ef107cf5SPeter Brune if (snes->nwork < 3) { 41*ef107cf5SPeter Brune ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 42*ef107cf5SPeter Brune } 43*ef107cf5SPeter Brune W = snes->work[0]; 44*ef107cf5SPeter Brune G = snes->work[1]; 45*ef107cf5SPeter Brune F = snes->work[2]; 46*ef107cf5SPeter Brune ierr = SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 47*ef107cf5SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 48*ef107cf5SPeter Brune ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr); 49*ef107cf5SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr); 50*ef107cf5SPeter Brune if (!colorcontainer) { 51*ef107cf5SPeter Brune /* create the coloring */ 52*ef107cf5SPeter Brune ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr); 53*ef107cf5SPeter Brune if (flg) { 54*ef107cf5SPeter Brune ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr); 55*ef107cf5SPeter Brune } else { 56*ef107cf5SPeter Brune if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);} 57*ef107cf5SPeter Brune ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr); 58*ef107cf5SPeter Brune ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr); 59*ef107cf5SPeter Brune ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 60*ef107cf5SPeter Brune ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr); 61*ef107cf5SPeter Brune ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 62*ef107cf5SPeter Brune } 63*ef107cf5SPeter Brune ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr); 64*ef107cf5SPeter Brune ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr); 65*ef107cf5SPeter Brune ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))GSDestroy_Private);CHKERRQ(ierr); 66*ef107cf5SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr); 67*ef107cf5SPeter Brune ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr); 68*ef107cf5SPeter Brune } else { 69*ef107cf5SPeter Brune ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr); 70*ef107cf5SPeter Brune } 71*ef107cf5SPeter Brune ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr); 72*ef107cf5SPeter Brune for (i=0;i<ncolors;i++) { 73*ef107cf5SPeter Brune for (k=0;k<its;k++) { 74*ef107cf5SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 75*ef107cf5SPeter Brune if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 76*ef107cf5SPeter Brune ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr); 77*ef107cf5SPeter Brune ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr); 78*ef107cf5SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 79*ef107cf5SPeter Brune for (j=0;j<size;j++) { 80*ef107cf5SPeter Brune ierr = VecSetValue(W,idx[j],h,ADD_VALUES);CHKERRQ(ierr); 81*ef107cf5SPeter Brune } 82*ef107cf5SPeter Brune ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr); 83*ef107cf5SPeter Brune if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);} 84*ef107cf5SPeter Brune for (j=0;j<size;j++) { 85*ef107cf5SPeter Brune ierr = VecGetValues(F,1,&idx[j],&f);CHKERRQ(ierr); 86*ef107cf5SPeter Brune ierr = VecGetValues(X,1,&idx[j],&x);CHKERRQ(ierr); 87*ef107cf5SPeter Brune ierr = VecGetValues(G,1,&idx[j],&g);CHKERRQ(ierr); 88*ef107cf5SPeter Brune ierr = VecGetValues(W,1,&idx[j],&w);CHKERRQ(ierr); 89*ef107cf5SPeter Brune d = (x*g-w*f) / PetscRealPart(g-f); 90*ef107cf5SPeter Brune ierr = VecSetValue(X,idx[j],d,INSERT_VALUES);CHKERRQ(ierr); 91*ef107cf5SPeter Brune } 92*ef107cf5SPeter Brune } 93*ef107cf5SPeter Brune } 94*ef107cf5SPeter Brune ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr); 95*ef107cf5SPeter Brune PetscFunctionReturn(0); 96*ef107cf5SPeter Brune } 97