Actual source code: gssecant.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <../src/snes/impls/gs/gsimpl.h>

  5: PetscErrorCode GSDestroy_Private(ISColoring coloring)
  6: {

 10:   ISColoringDestroy(&coloring);
 11:   return(0);

 13: }

 17: PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES snes,Vec X,Vec B,void *ctx)
 18: {
 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=0;
 31:   const PetscInt *idx;
 32:   PetscInt       size,s;
 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,isdone,alldone;
 39:   PetscScalar    *xa,*fa,*wa,*ga;

 42:   if (snes->nwork < 3) {
 43:     SNESSetWorkVecs(snes,3);
 44:   }
 45:   W = snes->work[0];
 46:   G = snes->work[1];
 47:   F = snes->work[2];
 48:   VecGetOwnershipRange(X,&s,NULL);
 49:   SNESGSGetTolerances(snes,&atol,&rtol,&stol,&its);
 50:   SNESGetDM(snes,&dm);
 51:   SNESGetFunction(snes,NULL,&func,&fctx);
 52:   PetscObjectQuery((PetscObject)snes,"SNESGSColoring",(PetscObject*)&colorcontainer);
 53:   if (!colorcontainer) {
 54:     /* create the coloring */
 55:     DMHasColoring(dm,&flg);
 56:     if (flg && !mat) {
 57:       DMCreateColoring(dm,IS_COLORING_GLOBAL,&coloring);
 58:     } else {
 59:       if (!snes->jacobian) {SNESSetUpMatrices(snes);}
 60:       MatColoringCreate(snes->jacobian,&mc);
 61:       MatColoringSetDistance(mc,1);
 62:       MatColoringSetFromOptions(mc);
 63:       MatColoringApply(mc,&coloring);
 64:       MatColoringDestroy(&mc);
 65:     }
 66:     PetscContainerCreate(PetscObjectComm((PetscObject)snes),&colorcontainer);
 67:     PetscContainerSetPointer(colorcontainer,(void *)coloring);
 68:     PetscContainerSetUserDestroy(colorcontainer,(PetscErrorCode (*)(void *))GSDestroy_Private);
 69:     PetscObjectCompose((PetscObject)snes,"SNESGSColoring",(PetscObject)colorcontainer);
 70:     PetscContainerDestroy(&colorcontainer);
 71:   } else {
 72:     PetscContainerGetPointer(colorcontainer,(void **)&coloring);
 73:   }
 74:   ISColoringGetIS(coloring,&ncolors,&coloris);
 75:   VecEqual(X,snes->vec_sol,&equal);
 76:   if (equal && snes->normschedule == SNES_NORM_ALWAYS) {
 77:     /* assume that the function is already computed */
 78:     VecCopy(snes->vec_func,F);
 79:   } else {
 80:     PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);
 81:     (*func)(snes,X,F,fctx);
 82:     PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);
 83:     if (B) {VecAXPY(F,-1.0,B);}
 84:   }
 85:   VecGetArray(X,&xa);
 86:   VecGetArray(F,&fa);
 87:   VecGetArray(G,&ga);
 88:   VecGetArray(W,&wa);
 89:   for (i=0;i<ncolors;i++) {
 90:     ISGetIndices(coloris[i],&idx);
 91:     ISGetLocalSize(coloris[i],&size);
 92:     VecCopy(X,W);
 93:     for (j=0;j<size;j++) {
 94:       wa[idx[j]-s] += h;
 95:     }
 96:     PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);
 97:     (*func)(snes,W,G,fctx);
 98:     PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);
 99:     if (B) {VecAXPY(G,-1.0,B);}
100:     for (k=0;k<its;k++) {
101:       dxt = 0.;
102:       xt = 0.;
103:       ft = 0.;
104:       for (j=0;j<size;j++) {
105:         f = fa[idx[j]-s];
106:         x = xa[idx[j]-s];
107:         g = ga[idx[j]-s];
108:         w = wa[idx[j]-s];
109:         if (PetscAbsScalar(g-f) > atol) {
110:           d = (x*g-w*f) / PetscRealPart(g-f);
111:         } else {
112:           d = x;
113:         }
114:         dxt += PetscRealPart(PetscSqr(d-x));
115:         xt += PetscRealPart(PetscSqr(x));
116:         ft += PetscRealPart(PetscSqr(f));
117:         xa[idx[j]-s] = d;
118:       }

120:       if (k == 0) ft1 = PetscSqrtReal(ft);
121:       if (k<its-1) {
122:         isdone = PETSC_FALSE;
123:         if (stol*PetscSqrtReal(xt) > PetscSqrtReal(dxt)) isdone = PETSC_TRUE;
124:         if (PetscSqrtReal(ft) < atol) isdone = PETSC_TRUE;
125:         if (rtol*ft1 > PetscSqrtReal(ft)) isdone = PETSC_TRUE;
126:         MPI_Allreduce(&isdone,&alldone,1,MPIU_BOOL,MPI_BAND,PetscObjectComm((PetscObject)snes));
127:         if (alldone) break;
128:       }
129:       if (i < ncolors-1 || k < its-1) {
130:         PetscLogEventBegin(SNES_GSFuncEval,snes,X,B,0);
131:         (*func)(snes,X,F,fctx);
132:         PetscLogEventEnd(SNES_GSFuncEval,snes,X,B,0);
133:         if (B) {VecAXPY(F,-1.0,B);}
134:       }
135:       if (k<its-1) {
136:         VecSwap(X,W);
137:         VecSwap(F,G);
138:       }
139:     }
140:   }
141:   VecRestoreArray(X,&xa);
142:   VecRestoreArray(F,&fa);
143:   VecRestoreArray(G,&ga);
144:   VecRestoreArray(W,&wa);
145:   ISColoringRestoreIS(coloring,&coloris);
146:   return(0);
147: }