1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> 2ef107cf5SPeter Brune 3ef107cf5SPeter Brune #undef __FUNCT__ 4ef107cf5SPeter Brune #define __FUNCT__ "GSDestroy_Private" 5ef107cf5SPeter Brune PetscErrorCode GSDestroy_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 15ef107cf5SPeter Brune #undef __FUNCT__ 16ef107cf5SPeter Brune #define __FUNCT__ "SNESComputeGSDefaultSecant" 17ef107cf5SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx) 18ef107cf5SPeter Brune { 19ef107cf5SPeter Brune PetscErrorCode ierr; 20737a7e12SPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 21ef107cf5SPeter Brune PetscInt i,j,k,ncolors; 22ef107cf5SPeter Brune DM dm; 23ef107cf5SPeter Brune PetscBool flg; 24ef107cf5SPeter Brune ISColoring coloring; 25ef107cf5SPeter Brune MatColoring mc; 26ef107cf5SPeter Brune Vec W,G,F; 27737a7e12SPeter Brune PetscScalar h=gs->h; 28ef107cf5SPeter Brune IS *coloris; 29ef107cf5SPeter Brune PetscScalar f,g,x,w,d; 30737a7e12SPeter Brune PetscReal dxt,xt,ft,ft1; 31ef107cf5SPeter Brune const PetscInt *idx; 32c03df580SPeter Brune PetscInt size,s; 33ef107cf5SPeter Brune PetscReal atol,rtol,stol; 34ef107cf5SPeter Brune PetscInt its; 35ef107cf5SPeter Brune PetscErrorCode (*func)(SNES,Vec,Vec,void*); 36ef107cf5SPeter Brune void *fctx; 37ef107cf5SPeter Brune PetscContainer colorcontainer; 38c03df580SPeter Brune PetscBool mat = gs->secant_mat,equal,isdone,alldone; 39c03df580SPeter Brune PetscScalar *xa,*fa,*wa,*ga; 40ef107cf5SPeter Brune 41ef107cf5SPeter Brune PetscFunctionBegin; 42ef107cf5SPeter Brune if (snes->nwork < 3) { 43ef107cf5SPeter Brune ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 44ef107cf5SPeter Brune } 45ef107cf5SPeter Brune W = snes->work[0]; 46ef107cf5SPeter Brune G = snes->work[1]; 47ef107cf5SPeter Brune F = snes->work[2]; 48c03df580SPeter Brune ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr); 49ef107cf5SPeter Brune ierr = SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 50ef107cf5SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 51ef107cf5SPeter Brune ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr); 52ef107cf5SPeter Brune ierr = PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr); 53ef107cf5SPeter Brune if (!colorcontainer) { 54ef107cf5SPeter Brune /* create the coloring */ 55ef107cf5SPeter Brune ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr); 56737a7e12SPeter Brune if (flg && !mat) { 57ef107cf5SPeter Brune ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr); 58ef107cf5SPeter Brune } else { 59ef107cf5SPeter Brune if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);} 60ef107cf5SPeter Brune ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr); 61ef107cf5SPeter Brune ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr); 62ef107cf5SPeter Brune ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 63ef107cf5SPeter Brune ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr); 64ef107cf5SPeter Brune ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 65ef107cf5SPeter Brune } 66ef107cf5SPeter Brune ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr); 67ef107cf5SPeter Brune ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr); 68ef107cf5SPeter Brune ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))GSDestroy_Private);CHKERRQ(ierr); 69ef107cf5SPeter Brune ierr = PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr); 70ef107cf5SPeter Brune ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr); 71ef107cf5SPeter Brune } else { 72ef107cf5SPeter Brune ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr); 73ef107cf5SPeter Brune } 74ef107cf5SPeter Brune ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr); 75737a7e12SPeter Brune ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr); 76c03df580SPeter Brune if (equal && snes->normschedule == SNES_NORM_ALWAYS) { 77737a7e12SPeter Brune /* assume that the function is already computed */ 78737a7e12SPeter Brune ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr); 79737a7e12SPeter Brune } else { 80e135a767SPeter Brune ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 81ef107cf5SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 82e135a767SPeter Brune ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 83ef107cf5SPeter Brune if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);} 84737a7e12SPeter Brune } 85c03df580SPeter Brune ierr = VecGetArray(X,&xa);CHKERRQ(ierr); 86c03df580SPeter Brune ierr = VecGetArray(F,&fa);CHKERRQ(ierr); 87c03df580SPeter Brune ierr = VecGetArray(G,&ga);CHKERRQ(ierr); 88c03df580SPeter Brune ierr = VecGetArray(W,&wa);CHKERRQ(ierr); 89737a7e12SPeter Brune for (i=0;i<ncolors;i++) { 90ef107cf5SPeter Brune ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr); 91ef107cf5SPeter Brune ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr); 92ef107cf5SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 93ef107cf5SPeter Brune for (j=0;j<size;j++) { 94c03df580SPeter Brune wa[idx[j]-s] += h; 95ef107cf5SPeter Brune } 96e135a767SPeter Brune ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 97ef107cf5SPeter Brune ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr); 98e135a767SPeter Brune ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 99ef107cf5SPeter Brune if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);} 100737a7e12SPeter Brune for (k=0;k<its;k++) { 101737a7e12SPeter Brune dxt = 0.; 102737a7e12SPeter Brune xt = 0.; 103c03df580SPeter Brune ft = 0.; 104ef107cf5SPeter Brune for (j=0;j<size;j++) { 105c03df580SPeter Brune f = fa[idx[j]-s]; 106c03df580SPeter Brune x = xa[idx[j]-s]; 107c03df580SPeter Brune g = ga[idx[j]-s]; 108c03df580SPeter Brune w = wa[idx[j]-s]; 109737a7e12SPeter Brune if (PetscAbsScalar(g-f) > atol) { 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 120*315bab25SPeter 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) { 130e135a767SPeter Brune ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr); 131737a7e12SPeter Brune ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr); 132e135a767SPeter Brune ierr = PetscLogEventEnd(SNES_GSFuncEval,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