xref: /petsc/src/snes/impls/gs/gssecant.c (revision 737a7e128d3ed57ea128fd0d033fca43bac7efb0)
1*737a7e12SPeter 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;
20*737a7e12SPeter 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;
27*737a7e12SPeter Brune   PetscScalar    h=gs->h;
28ef107cf5SPeter Brune   IS             *coloris;
29ef107cf5SPeter Brune   PetscScalar    f,g,x,w,d;
30*737a7e12SPeter Brune   PetscReal      dxt,xt,ft,ft1;
31ef107cf5SPeter Brune   const PetscInt *idx;
32ef107cf5SPeter Brune   PetscInt       size;
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;
38*737a7e12SPeter Brune   PetscBool      mat = gs->secant_mat,equal;
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];
47ef107cf5SPeter Brune   ierr = SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr);
48ef107cf5SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
49ef107cf5SPeter Brune   ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr);
50ef107cf5SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr);
51ef107cf5SPeter Brune   if (!colorcontainer) {
52ef107cf5SPeter Brune     /* create the coloring */
53ef107cf5SPeter Brune     ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr);
54*737a7e12SPeter Brune     if (flg && !mat) {
55ef107cf5SPeter Brune       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr);
56ef107cf5SPeter Brune     } else {
57ef107cf5SPeter Brune       if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
58ef107cf5SPeter Brune       ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr);
59ef107cf5SPeter Brune       ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
60ef107cf5SPeter Brune       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
61ef107cf5SPeter Brune       ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr);
62ef107cf5SPeter Brune       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
63ef107cf5SPeter Brune     }
64ef107cf5SPeter Brune     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr);
65ef107cf5SPeter Brune     ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr);
66ef107cf5SPeter Brune     ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))GSDestroy_Private);CHKERRQ(ierr);
67ef107cf5SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr);
68ef107cf5SPeter Brune     ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr);
69ef107cf5SPeter Brune   } else {
70ef107cf5SPeter Brune     ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr);
71ef107cf5SPeter Brune   }
72ef107cf5SPeter Brune   ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr);
73*737a7e12SPeter Brune   ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr);
74*737a7e12SPeter Brune   if (equal) {
75*737a7e12SPeter Brune     /* assume that the function is already computed */
76*737a7e12SPeter Brune     ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr);
77*737a7e12SPeter Brune   } else {
78ef107cf5SPeter Brune     ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
79ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
80*737a7e12SPeter Brune   }
81*737a7e12SPeter Brune   for (i=0;i<ncolors;i++) {
82ef107cf5SPeter Brune     ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr);
83ef107cf5SPeter Brune     ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr);
84ef107cf5SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
85ef107cf5SPeter Brune     for (j=0;j<size;j++) {
86ef107cf5SPeter Brune       ierr = VecSetValue(W,idx[j],h,ADD_VALUES);CHKERRQ(ierr);
87ef107cf5SPeter Brune     }
88ef107cf5SPeter Brune     ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr);
89ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);}
90*737a7e12SPeter Brune     for (k=0;k<its;k++) {
91*737a7e12SPeter Brune       dxt = 0.;
92*737a7e12SPeter Brune       xt = 0.;
93ef107cf5SPeter Brune       for (j=0;j<size;j++) {
94ef107cf5SPeter Brune         ierr = VecGetValues(F,1,&idx[j],&f);CHKERRQ(ierr);
95ef107cf5SPeter Brune         ierr = VecGetValues(X,1,&idx[j],&x);CHKERRQ(ierr);
96ef107cf5SPeter Brune         ierr = VecGetValues(G,1,&idx[j],&g);CHKERRQ(ierr);
97ef107cf5SPeter Brune         ierr = VecGetValues(W,1,&idx[j],&w);CHKERRQ(ierr);
98*737a7e12SPeter Brune         if (PetscAbsScalar(g-f) > atol) {
99ef107cf5SPeter Brune           d = (x*g-w*f) / PetscRealPart(g-f);
100*737a7e12SPeter Brune         } else {
101*737a7e12SPeter Brune           d = x;
102*737a7e12SPeter Brune         }
103*737a7e12SPeter Brune         dxt += (d-x)*(d-x);
104*737a7e12SPeter Brune         xt += x*x;
105*737a7e12SPeter Brune         ft += f*f;
106ef107cf5SPeter Brune         ierr = VecSetValue(X,idx[j],d,INSERT_VALUES);CHKERRQ(ierr);
107ef107cf5SPeter Brune       }
108*737a7e12SPeter Brune       if (k == 0) ft1 = PetscSqrtScalar(ft);
109*737a7e12SPeter Brune       if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) break;
110*737a7e12SPeter Brune       if (PetscSqrtReal(ft) < atol) break;
111*737a7e12SPeter Brune       if (rtol*ft1 > PetscSqrtReal(ft)) break;
112*737a7e12SPeter Brune       if (i < ncolors-1 || k < its-1) {
113*737a7e12SPeter Brune         ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
114*737a7e12SPeter Brune         if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
115*737a7e12SPeter Brune       }
116*737a7e12SPeter Brune       if (k<its-1) {
117*737a7e12SPeter Brune         ierr = VecSwap(X,W);CHKERRQ(ierr);
118*737a7e12SPeter Brune         ierr = VecSwap(F,G);CHKERRQ(ierr);
119*737a7e12SPeter Brune       }
120ef107cf5SPeter Brune     }
121ef107cf5SPeter Brune   }
122ef107cf5SPeter Brune   ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr);
123ef107cf5SPeter Brune   PetscFunctionReturn(0);
124ef107cf5SPeter Brune }
125