xref: /petsc/src/snes/impls/gs/gssecant.c (revision 235fd6e69a2bea97b8017ef133f6bee539b5dc6b)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h>
2ef107cf5SPeter Brune 
3ef107cf5SPeter Brune #undef __FUNCT__
4*235fd6e6SBarry Smith #define __FUNCT__ "SNESNGSDestroy_Private"
5*235fd6e6SBarry Smith static PetscErrorCode SNESNGSDestroy_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 #undef __FUNCT__
15ef107cf5SPeter Brune #define __FUNCT__ "SNESComputeGSDefaultSecant"
16ef107cf5SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx)
17ef107cf5SPeter Brune {
18ef107cf5SPeter Brune   PetscErrorCode ierr;
19737a7e12SPeter Brune   SNES_GS        *gs = (SNES_GS*)snes->data;
20ef107cf5SPeter Brune   PetscInt       i,j,k,ncolors;
21ef107cf5SPeter Brune   DM             dm;
22ef107cf5SPeter Brune   PetscBool      flg;
23ef107cf5SPeter Brune   ISColoring     coloring;
24ef107cf5SPeter Brune   MatColoring    mc;
25ef107cf5SPeter Brune   Vec            W,G,F;
26737a7e12SPeter Brune   PetscScalar    h=gs->h;
27ef107cf5SPeter Brune   IS             *coloris;
28ef107cf5SPeter Brune   PetscScalar    f,g,x,w,d;
29c1132020SPeter Brune   PetscReal      dxt,xt,ft,ft1=0;
30ef107cf5SPeter Brune   const PetscInt *idx;
31c03df580SPeter Brune   PetscInt       size,s;
32ef107cf5SPeter Brune   PetscReal      atol,rtol,stol;
33ef107cf5SPeter Brune   PetscInt       its;
34ef107cf5SPeter Brune   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
35ef107cf5SPeter Brune   void           *fctx;
36ef107cf5SPeter Brune   PetscContainer colorcontainer;
37c03df580SPeter Brune   PetscBool      mat = gs->secant_mat,equal,isdone,alldone;
38c03df580SPeter Brune   PetscScalar    *xa,*fa,*wa,*ga;
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];
47c03df580SPeter Brune   ierr = VecGetOwnershipRange(X,&s,NULL);CHKERRQ(ierr);
48ef107cf5SPeter Brune   ierr = SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr);
49ef107cf5SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
50ef107cf5SPeter Brune   ierr = SNESGetFunction(snes,NULL,&func,&fctx);CHKERRQ(ierr);
51ef107cf5SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);CHKERRQ(ierr);
52ef107cf5SPeter Brune   if (!colorcontainer) {
53ef107cf5SPeter Brune     /* create the coloring */
54ef107cf5SPeter Brune     ierr = DMHasColoring(dm,&flg);CHKERRQ(ierr);
55737a7e12SPeter Brune     if (flg && !mat) {
56ef107cf5SPeter Brune       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);CHKERRQ(ierr);
57ef107cf5SPeter Brune     } else {
58ef107cf5SPeter Brune       if (!snes->jacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
59ef107cf5SPeter Brune       ierr = MatColoringCreate(snes->jacobian,&mc);CHKERRQ(ierr);
60ef107cf5SPeter Brune       ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
61ef107cf5SPeter Brune       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
62ef107cf5SPeter Brune       ierr = MatColoringApply(mc,&coloring);CHKERRQ(ierr);
63ef107cf5SPeter Brune       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
64ef107cf5SPeter Brune     }
65ef107cf5SPeter Brune     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);CHKERRQ(ierr);
66ef107cf5SPeter Brune     ierr = PetscContainerSetPointer(colorcontainer,(void *)coloring);CHKERRQ(ierr);
67*235fd6e6SBarry Smith     ierr = PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))SNESNGSDestroy_Private);CHKERRQ(ierr);
68ef107cf5SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);CHKERRQ(ierr);
69ef107cf5SPeter Brune     ierr = PetscContainerDestroy(&colorcontainer);CHKERRQ(ierr);
70ef107cf5SPeter Brune   } else {
71ef107cf5SPeter Brune     ierr = PetscContainerGetPointer(colorcontainer,(void **)&coloring);CHKERRQ(ierr);
72ef107cf5SPeter Brune   }
73ef107cf5SPeter Brune   ierr = ISColoringGetIS(coloring,&ncolors,&coloris);CHKERRQ(ierr);
74737a7e12SPeter Brune   ierr = VecEqual(X,snes->vec_sol,&equal);CHKERRQ(ierr);
75c03df580SPeter Brune   if (equal && snes->normschedule == SNES_NORM_ALWAYS) {
76737a7e12SPeter Brune     /* assume that the function is already computed */
77737a7e12SPeter Brune     ierr = VecCopy(snes->vec_func,F);CHKERRQ(ierr);
78737a7e12SPeter Brune   } else {
79e135a767SPeter Brune     ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr);
80ef107cf5SPeter Brune     ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
81e135a767SPeter Brune     ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr);
82ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
83737a7e12SPeter Brune   }
84c03df580SPeter Brune   ierr = VecGetArray(X,&xa);CHKERRQ(ierr);
85c03df580SPeter Brune   ierr = VecGetArray(F,&fa);CHKERRQ(ierr);
86c03df580SPeter Brune   ierr = VecGetArray(G,&ga);CHKERRQ(ierr);
87c03df580SPeter Brune   ierr = VecGetArray(W,&wa);CHKERRQ(ierr);
88737a7e12SPeter Brune   for (i=0;i<ncolors;i++) {
89ef107cf5SPeter Brune     ierr = ISGetIndices(coloris[i],&idx);CHKERRQ(ierr);
90ef107cf5SPeter Brune     ierr = ISGetLocalSize(coloris[i],&size);CHKERRQ(ierr);
91ef107cf5SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
92ef107cf5SPeter Brune     for (j=0;j<size;j++) {
93c03df580SPeter Brune       wa[idx[j]-s] += h;
94ef107cf5SPeter Brune     }
95e135a767SPeter Brune     ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr);
96ef107cf5SPeter Brune     ierr = (*func)(snes,W,G,fctx);CHKERRQ(ierr);
97e135a767SPeter Brune     ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr);
98ef107cf5SPeter Brune     if (B) {ierr = VecAXPY(G,-1.0,B);CHKERRQ(ierr);}
99737a7e12SPeter Brune     for (k=0;k<its;k++) {
100737a7e12SPeter Brune       dxt = 0.;
101737a7e12SPeter Brune       xt = 0.;
102c03df580SPeter Brune       ft = 0.;
103ef107cf5SPeter Brune       for (j=0;j<size;j++) {
104c03df580SPeter Brune         f = fa[idx[j]-s];
105c03df580SPeter Brune         x = xa[idx[j]-s];
106c03df580SPeter Brune         g = ga[idx[j]-s];
107c03df580SPeter Brune         w = wa[idx[j]-s];
108737a7e12SPeter Brune         if (PetscAbsScalar(g-f) > atol) {
109ef107cf5SPeter Brune           d = (x*g-w*f) / PetscRealPart(g-f);
110737a7e12SPeter Brune         } else {
111737a7e12SPeter Brune           d = x;
112737a7e12SPeter Brune         }
113400ea6b8SPeter Brune         dxt += PetscRealPart(PetscSqr(d-x));
114400ea6b8SPeter Brune         xt += PetscRealPart(PetscSqr(x));
115400ea6b8SPeter Brune         ft += PetscRealPart(PetscSqr(f));
116c03df580SPeter Brune         xa[idx[j]-s] = d;
117ef107cf5SPeter Brune       }
118c03df580SPeter Brune 
119315bab25SPeter Brune       if (k == 0) ft1 = PetscSqrtReal(ft);
120c03df580SPeter Brune       if (k<its-1) {
121c03df580SPeter Brune         isdone = PETSC_FALSE;
122c03df580SPeter Brune         if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE;
123c03df580SPeter Brune         if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE;
124c03df580SPeter Brune         if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE;
125c03df580SPeter Brune         ierr = MPI_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
126c03df580SPeter Brune         if (alldone) break;
127c03df580SPeter Brune       }
128737a7e12SPeter Brune       if (i < ncolors-1 || k < its-1) {
129e135a767SPeter Brune         ierr = PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr);
130737a7e12SPeter Brune         ierr = (*func)(snes,X,F,fctx);CHKERRQ(ierr);
131e135a767SPeter Brune         ierr = PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);CHKERRQ(ierr);
132737a7e12SPeter Brune         if (B) {ierr = VecAXPY(F,-1.0,B);CHKERRQ(ierr);}
133737a7e12SPeter Brune       }
134737a7e12SPeter Brune       if (k<its-1) {
135737a7e12SPeter Brune         ierr = VecSwap(X,W);CHKERRQ(ierr);
136737a7e12SPeter Brune         ierr = VecSwap(F,G);CHKERRQ(ierr);
137737a7e12SPeter Brune       }
138ef107cf5SPeter Brune     }
139ef107cf5SPeter Brune   }
140c03df580SPeter Brune   ierr = VecRestoreArray(X,&xa);CHKERRQ(ierr);
141c03df580SPeter Brune   ierr = VecRestoreArray(F,&fa);CHKERRQ(ierr);
142c03df580SPeter Brune   ierr = VecRestoreArray(G,&ga);CHKERRQ(ierr);
143c03df580SPeter Brune   ierr = VecRestoreArray(W,&wa);CHKERRQ(ierr);
144ef107cf5SPeter Brune   ierr = ISColoringRestoreIS(coloring,&coloris);CHKERRQ(ierr);
145ef107cf5SPeter Brune   PetscFunctionReturn(0);
146ef107cf5SPeter Brune }
147