xref: /petsc/src/snes/impls/gs/gssecant.c (revision 071fcb05f2a6aea8aef2c2090580530b9e9df77e)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h>
2ef107cf5SPeter Brune 
3be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx)
4ef107cf5SPeter Brune {
5ef107cf5SPeter Brune   PetscErrorCode ierr;
6be95d8f1SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
7ef107cf5SPeter Brune   PetscInt       i,j,k,ncolors;
8ef107cf5SPeter Brune   DM             dm;
9ef107cf5SPeter Brune   PetscBool      flg;
108a86d3c5SBarry Smith   ISColoring     coloring = gs->coloring;
11ef107cf5SPeter Brune   MatColoring    mc;
12ef107cf5SPeter Brune   Vec            W,G,F;
13737a7e12SPeter Brune   PetscScalar    h=gs->h;
14ef107cf5SPeter Brune   IS             *coloris;
15ef107cf5SPeter Brune   PetscScalar    f,g,x,w,d;
16c1132020SPeter Brune   PetscReal      dxt,xt,ft,ft1=0;
17ef107cf5SPeter Brune   const PetscInt *idx;
18c03df580SPeter Brune   PetscInt       size,s;
19ef107cf5SPeter Brune   PetscReal      atol,rtol,stol;
20ef107cf5SPeter Brune   PetscInt       its;
21ef107cf5SPeter Brune   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
22ef107cf5SPeter Brune   void           *fctx;
23c03df580SPeter Brune   PetscBool      mat = gs->secant_mat,equal,isdone,alldone;
24c03df580SPeter Brune   PetscScalar    *xa,*fa,*wa,*ga;
25ef107cf5SPeter Brune 
26ef107cf5SPeter Brune   PetscFunctionBegin;
27ef107cf5SPeter Brune   if (snes->nwork < 3) {
28ef107cf5SPeter Brune     ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
29ef107cf5SPeter Brune   }
30ef107cf5SPeter Brune   W = snes->work[0];
31ef107cf5SPeter Brune   G = snes->work[1];
32ef107cf5SPeter Brune   F = snes->work[2];
33c03df580SPeter Brune   ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr);
34be95d8f1SBarry Smith   ierr = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr);
35ef107cf5SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
36ef107cf5SPeter Brune   ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr);
378a86d3c5SBarry Smith   if (!coloring) {
38ef107cf5SPeter Brune     /* create the coloring */
39ef107cf5SPeter Brune     ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr);
40737a7e12SPeter Brune     if (flg && !mat) {
41ef107cf5SPeter Brune       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr);
42ef107cf5SPeter Brune     } else {
43ef107cf5SPeter Brune       if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
44ef107cf5SPeter Brune       ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr);
45ef107cf5SPeter Brune       ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
46ef107cf5SPeter Brune       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
47ef107cf5SPeter Brune       ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr);
48ef107cf5SPeter Brune       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
49ef107cf5SPeter Brune     }
508a86d3c5SBarry Smith     gs->coloring = coloring;
51ef107cf5SPeter Brune   }
52*071fcb05SBarry Smith   ierr = ISColoringGetIS(coloring,PETSC_USE_POINTER,&ncolors,&coloris);CHKERRQ(ierr);
53737a7e12SPeter Brune   ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr);
54c03df580SPeter Brune   if (equal && snes->normschedule == SNES_NORM_ALWAYS) {
55737a7e12SPeter Brune     /* assume that the function is already computed */
56737a7e12SPeter Brune     ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr);
57737a7e12SPeter Brune   } else {
58be95d8f1SBarry Smith     ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
59ef107cf5SPeter Brune     ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
60be95d8f1SBarry Smith     ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
61ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
62737a7e12SPeter Brune   }
63c03df580SPeter Brune   ierr = VecGetArray(X,&xa);CHKERRQ(ierr);
64c03df580SPeter Brune   ierr = VecGetArray(F,&fa);CHKERRQ(ierr);
65c03df580SPeter Brune   ierr = VecGetArray(G,&ga);CHKERRQ(ierr);
66c03df580SPeter Brune   ierr = VecGetArray(W,&wa);CHKERRQ(ierr);
67737a7e12SPeter Brune   for (i=0;i<ncolors;i++) {
68ef107cf5SPeter Brune     ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr);
69ef107cf5SPeter Brune     ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr);
70ef107cf5SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
71ef107cf5SPeter Brune     for (j=0;j<size;j++) {
72c03df580SPeter Brune       wa[idx[j]-s] += h;
73ef107cf5SPeter Brune     }
74be95d8f1SBarry Smith     ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
75ef107cf5SPeter Brune     ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr);
76be95d8f1SBarry Smith     ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
77ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);}
78737a7e12SPeter Brune     for (k=0;k<its;k++) {
79737a7e12SPeter Brune       dxt = 0.;
80737a7e12SPeter Brune       xt = 0.;
81c03df580SPeter Brune       ft = 0.;
82ef107cf5SPeter Brune       for (j=0;j<size;j++) {
83c03df580SPeter Brune         f = fa[idx[j]-s];
84c03df580SPeter Brune         x = xa[idx[j]-s];
85c03df580SPeter Brune         g = ga[idx[j]-s];
86c03df580SPeter Brune         w = wa[idx[j]-s];
87737a7e12SPeter Brune         if (PetscAbsScalar(g-f) > atol) {
884b655d89SMatthew G. Knepley           /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */
89ef107cf5SPeter Brune           d = (x*g-w*f) / PetscRealPart(g-f);
90737a7e12SPeter Brune         } else {
91737a7e12SPeter Brune           d = x;
92737a7e12SPeter Brune         }
93400ea6b8SPeter Brune         dxt += PetscRealPart(PetscSqr(d-x));
94400ea6b8SPeter Brune         xt += PetscRealPart(PetscSqr(x));
95400ea6b8SPeter Brune         ft += PetscRealPart(PetscSqr(f));
96c03df580SPeter Brune         xa[idx[j]-s] = d;
97ef107cf5SPeter Brune       }
98c03df580SPeter Brune 
99315bab25SPeter Brune       if (k == 0) ft1 = PetscSqrtReal(ft);
100c03df580SPeter Brune       if (k<its-1) {
101c03df580SPeter Brune         isdone = PETSC_FALSE;
102c03df580SPeter Brune         if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE;
103c03df580SPeter Brune         if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE;
104c03df580SPeter Brune         if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE;
105b2566f29SBarry Smith         ierr = MPIU_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
106c03df580SPeter Brune         if (alldone) break;
107c03df580SPeter Brune       }
108737a7e12SPeter Brune       if (i < ncolors-1 || k < its-1) {
109be95d8f1SBarry Smith         ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
110737a7e12SPeter Brune         ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
111be95d8f1SBarry Smith         ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
112737a7e12SPeter Brune         if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
113737a7e12SPeter Brune       }
114737a7e12SPeter Brune       if (k<its-1) {
115737a7e12SPeter Brune         ierr = VecSwap(X,W);CHKERRQ(ierr);
116737a7e12SPeter Brune         ierr = VecSwap(F,G);CHKERRQ(ierr);
117737a7e12SPeter Brune       }
118ef107cf5SPeter Brune     }
119ef107cf5SPeter Brune   }
120c03df580SPeter Brune   ierr = VecRestoreArray(X,&xa);CHKERRQ(ierr);
121c03df580SPeter Brune   ierr = VecRestoreArray(F,&fa);CHKERRQ(ierr);
122c03df580SPeter Brune   ierr = VecRestoreArray(G,&ga);CHKERRQ(ierr);
123c03df580SPeter Brune   ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr);
124*071fcb05SBarry Smith   ierr = ISColoringRestoreIS(coloring,PETSC_USE_POINTER,&coloris);CHKERRQ(ierr);
125ef107cf5SPeter Brune   PetscFunctionReturn(0);
126ef107cf5SPeter Brune }
127