xref: /petsc/src/snes/impls/gs/gssecant.c (revision 737a7e128d3ed57ea128fd0d033fca43bac7efb0)
1 #include <../src/snes/impls/gs/gsimpl.h>
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "GSDestroy_Private"
5 PetscErrorCode GSDestroy_Private(ISColoring coloring)
6 {
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
11   PetscFunctionReturn(0);
12 
13 }
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "SNESComputeGSDefaultSecant"
17 PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx)
18 {
19   PetscErrorCode ierr;
20   SNES_GS *gs = (SNES_GS*)snes->data;
21   PetscInt       i,j,k,ncolors;
22   DM             dm;
23   PetscBool      flg;
24   ISColoring     coloring;
25   MatColoring    mc;
26   Vec            W,G,F;
27   PetscScalar    h=gs->h;
28   IS             *coloris;
29   PetscScalar    f,g,x,w,d;
30   PetscReal      dxt,xt,ft,ft1;
31   const PetscInt *idx;
32   PetscInt       size;
33   PetscReal      atol,rtol,stol;
34   PetscInt       its;
35   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
36   void           *fctx;
37   PetscContainer colorcontainer;
38   PetscBool      mat = gs->secant_mat,equal;
39 
40   PetscFunctionBegin;
41   if (snes->nwork < 3) {
42     ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
43   }
44   W = snes->work[0];
45   G = snes->work[1];
46   F = snes->work[2];
47   ierr = SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr);
48   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
49   ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr);
50   ierr = PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr);
51   if (!colorcontainer) {
52     /* create the coloring */
53     ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr);
54     if (flg && !mat) {
55       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr);
56     } else {
57       if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
58       ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr);
59       ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
60       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
61       ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr);
62       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
63     }
64     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr);
65     ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr);
66     ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))GSDestroy_Private);CHKERRQ(ierr);
67     ierr = PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr);
68     ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr);
69   } else {
70     ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr);
71   }
72   ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr);
73   ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr);
74   if (equal) {
75     /* assume that the function is already computed */
76     ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr);
77   } else {
78     ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
79     if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
80   }
81   for (i=0;i<ncolors;i++) {
82     ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr);
83     ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr);
84     ierr = VecCopy(X,W);CHKERRQ(ierr);
85     for (j=0;j<size;j++) {
86       ierr = VecSetValue(W,idx[j],h,ADD_VALUES);CHKERRQ(ierr);
87     }
88     ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr);
89     if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);}
90     for (k=0;k<its;k++) {
91       dxt = 0.;
92       xt = 0.;
93       for (j=0;j<size;j++) {
94         ierr = VecGetValues(F,1,&idx[j],&f);CHKERRQ(ierr);
95         ierr = VecGetValues(X,1,&idx[j],&x);CHKERRQ(ierr);
96         ierr = VecGetValues(G,1,&idx[j],&g);CHKERRQ(ierr);
97         ierr = VecGetValues(W,1,&idx[j],&w);CHKERRQ(ierr);
98         if (PetscAbsScalar(g-f) > atol) {
99           d = (x*g-w*f) / PetscRealPart(g-f);
100         } else {
101           d = x;
102         }
103         dxt += (d-x)*(d-x);
104         xt += x*x;
105         ft += f*f;
106         ierr = VecSetValue(X,idx[j],d,INSERT_VALUES);CHKERRQ(ierr);
107       }
108       if (k == 0) ft1 = PetscSqrtScalar(ft);
109       if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) break;
110       if (PetscSqrtReal(ft) < atol) break;
111       if (rtol*ft1 > PetscSqrtReal(ft)) break;
112       if (i < ncolors-1 || k < its-1) {
113         ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
114         if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
115       }
116       if (k<its-1) {
117         ierr = VecSwap(X,W);CHKERRQ(ierr);
118         ierr = VecSwap(F,G);CHKERRQ(ierr);
119       }
120     }
121   }
122   ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr);
123   PetscFunctionReturn(0);
124 }
125