xref: /petsc/src/snes/impls/gs/gssecant.c (revision 8a86d3c58f8fa1b9ff9e69bfa0d33045bfc2518d)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h>
2ef107cf5SPeter Brune 
3ef107cf5SPeter Brune #undef __FUNCT__
4be95d8f1SBarry Smith #define __FUNCT__ "SNESComputeNGSDefaultSecant"
5be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx)
6ef107cf5SPeter Brune {
7ef107cf5SPeter Brune   PetscErrorCode ierr;
8be95d8f1SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
9ef107cf5SPeter Brune   PetscInt       i,j,k,ncolors;
10ef107cf5SPeter Brune   DM             dm;
11ef107cf5SPeter Brune   PetscBool      flg;
12*8a86d3c5SBarry Smith   ISColoring     coloring = gs->coloring;
13ef107cf5SPeter Brune   MatColoring    mc;
14ef107cf5SPeter Brune   Vec            W,G,F;
15737a7e12SPeter Brune   PetscScalar    h=gs->h;
16ef107cf5SPeter Brune   IS             *coloris;
17ef107cf5SPeter Brune   PetscScalar    f,g,x,w,d;
18c1132020SPeter Brune   PetscReal      dxt,xt,ft,ft1=0;
19ef107cf5SPeter Brune   const PetscInt *idx;
20c03df580SPeter Brune   PetscInt       size,s;
21ef107cf5SPeter Brune   PetscReal      atol,rtol,stol;
22ef107cf5SPeter Brune   PetscInt       its;
23ef107cf5SPeter Brune   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
24ef107cf5SPeter Brune   void           *fctx;
25c03df580SPeter Brune   PetscBool      mat = gs->secant_mat,equal,isdone,alldone;
26c03df580SPeter Brune   PetscScalar    *xa,*fa,*wa,*ga;
27ef107cf5SPeter Brune 
28ef107cf5SPeter Brune   PetscFunctionBegin;
29ef107cf5SPeter Brune   if (snes->nwork < 3) {
30ef107cf5SPeter Brune     ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
31ef107cf5SPeter Brune   }
32ef107cf5SPeter Brune   W = snes->work[0];
33ef107cf5SPeter Brune   G = snes->work[1];
34ef107cf5SPeter Brune   F = snes->work[2];
35c03df580SPeter Brune   ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr);
36be95d8f1SBarry Smith   ierr = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr);
37ef107cf5SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
38ef107cf5SPeter Brune   ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr);
39*8a86d3c5SBarry Smith   if (!coloring) {
40ef107cf5SPeter Brune     /* create the coloring */
41ef107cf5SPeter Brune     ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr);
42737a7e12SPeter Brune     if (flg && !mat) {
43ef107cf5SPeter Brune       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr);
44ef107cf5SPeter Brune     } else {
45ef107cf5SPeter Brune       if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
46ef107cf5SPeter Brune       ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr);
47ef107cf5SPeter Brune       ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
48ef107cf5SPeter Brune       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
49ef107cf5SPeter Brune       ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr);
50ef107cf5SPeter Brune       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
51ef107cf5SPeter Brune     }
52*8a86d3c5SBarry Smith     gs->coloring = coloring;
53ef107cf5SPeter Brune   }
54ef107cf5SPeter Brune   ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr);
55737a7e12SPeter Brune   ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr);
56c03df580SPeter Brune   if (equal && snes->normschedule == SNES_NORM_ALWAYS) {
57737a7e12SPeter Brune     /* assume that the function is already computed */
58737a7e12SPeter Brune     ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr);
59737a7e12SPeter Brune   } else {
60be95d8f1SBarry Smith     ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
61ef107cf5SPeter Brune     ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
62be95d8f1SBarry Smith     ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
63ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
64737a7e12SPeter Brune   }
65c03df580SPeter Brune   ierr = VecGetArray(X,&xa);CHKERRQ(ierr);
66c03df580SPeter Brune   ierr = VecGetArray(F,&fa);CHKERRQ(ierr);
67c03df580SPeter Brune   ierr = VecGetArray(G,&ga);CHKERRQ(ierr);
68c03df580SPeter Brune   ierr = VecGetArray(W,&wa);CHKERRQ(ierr);
69737a7e12SPeter Brune   for (i=0;i<ncolors;i++) {
70ef107cf5SPeter Brune     ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr);
71ef107cf5SPeter Brune     ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr);
72ef107cf5SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
73ef107cf5SPeter Brune     for (j=0;j<size;j++) {
74c03df580SPeter Brune       wa[idx[j]-s] += h;
75ef107cf5SPeter Brune     }
76be95d8f1SBarry Smith     ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
77ef107cf5SPeter Brune     ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr);
78be95d8f1SBarry Smith     ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
79ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);}
80737a7e12SPeter Brune     for (k=0;k<its;k++) {
81737a7e12SPeter Brune       dxt = 0.;
82737a7e12SPeter Brune       xt = 0.;
83c03df580SPeter Brune       ft = 0.;
84ef107cf5SPeter Brune       for (j=0;j<size;j++) {
85c03df580SPeter Brune         f = fa[idx[j]-s];
86c03df580SPeter Brune         x = xa[idx[j]-s];
87c03df580SPeter Brune         g = ga[idx[j]-s];
88c03df580SPeter Brune         w = wa[idx[j]-s];
89737a7e12SPeter Brune         if (PetscAbsScalar(g-f) > atol) {
904b655d89SMatthew G. Knepley           /* This is equivalent to d = x - (h*f) / PetscRealPart(g-f) */
91ef107cf5SPeter Brune           d = (x*g-w*f) / PetscRealPart(g-f);
92737a7e12SPeter Brune         } else {
93737a7e12SPeter Brune           d = x;
94737a7e12SPeter Brune         }
95400ea6b8SPeter Brune         dxt += PetscRealPart(PetscSqr(d-x));
96400ea6b8SPeter Brune         xt += PetscRealPart(PetscSqr(x));
97400ea6b8SPeter Brune         ft += PetscRealPart(PetscSqr(f));
98c03df580SPeter Brune         xa[idx[j]-s] = d;
99ef107cf5SPeter Brune       }
100c03df580SPeter Brune 
101315bab25SPeter Brune       if (k == 0) ft1 = PetscSqrtReal(ft);
102c03df580SPeter Brune       if (k<its-1) {
103c03df580SPeter Brune         isdone = PETSC_FALSE;
104c03df580SPeter Brune         if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE;
105c03df580SPeter Brune         if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE;
106c03df580SPeter Brune         if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE;
107b2566f29SBarry Smith         ierr = MPIU_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
108c03df580SPeter Brune         if (alldone) break;
109c03df580SPeter Brune       }
110737a7e12SPeter Brune       if (i < ncolors-1 || k < its-1) {
111be95d8f1SBarry Smith         ierr = PetscLogEventBegin(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
112737a7e12SPeter Brune         ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
113be95d8f1SBarry Smith         ierr = PetscLogEventEnd(SNES_NGSFuncEval,snes,X,B,0);CHKERRQ(ierr);
114737a7e12SPeter Brune         if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
115737a7e12SPeter Brune       }
116737a7e12SPeter Brune       if (k<its-1) {
117737a7e12SPeter Brune         ierr = VecSwap(X,W);CHKERRQ(ierr);
118737a7e12SPeter Brune         ierr = VecSwap(F,G);CHKERRQ(ierr);
119737a7e12SPeter Brune       }
120ef107cf5SPeter Brune     }
121ef107cf5SPeter Brune   }
122c03df580SPeter Brune   ierr = VecRestoreArray(X,&xa);CHKERRQ(ierr);
123c03df580SPeter Brune   ierr = VecRestoreArray(F,&fa);CHKERRQ(ierr);
124c03df580SPeter Brune   ierr = VecRestoreArray(G,&ga);CHKERRQ(ierr);
125c03df580SPeter Brune   ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr);
126ef107cf5SPeter Brune   ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr);
127ef107cf5SPeter Brune   PetscFunctionReturn(0);
128ef107cf5SPeter Brune }
129