Actual source code: snesnoise.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: #include <petsc/private/snesimpl.h>

  4: /* Data used by Jorge's diff parameter computation method */
  5: typedef struct {
  6:   Vec      *workv;           /* work vectors */
  7:   FILE     *fp;              /* output file */
  8:   int      function_count;   /* count of function evaluations for diff param estimation */
  9:   double   fnoise_min;       /* minimim allowable noise */
 10:   double   hopt_min;         /* minimum allowable hopt */
 11:   double   h_first_try;      /* first try for h used in diff parameter estimate */
 12:   PetscInt fnoise_resets;    /* number of times we've reset the noise estimate */
 13:   PetscInt hopt_resets;      /* number of times we've reset the hopt estimate */
 14: } DIFFPAR_MORE;


 17: extern PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*);
 18: static PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double);
 19: extern PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double);
 20: extern PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes);

 24: PetscErrorCode SNESDiffParameterCreate_More(SNES snes,Vec x,void **outneP)
 25: {
 26:   DIFFPAR_MORE   *neP;
 27:   Vec            w;
 28:   PetscRandom    rctx;  /* random number generator context */
 30:   PetscBool      flg;
 31:   char           noise_file[PETSC_MAX_PATH_LEN];

 34:   PetscNewLog(snes,&neP);

 36:   neP->function_count = 0;
 37:   neP->fnoise_min     = 1.0e-20;
 38:   neP->hopt_min       = 1.0e-8;
 39:   neP->h_first_try    = 1.0e-3;
 40:   neP->fnoise_resets  = 0;
 41:   neP->hopt_resets    = 0;

 43:   /* Create work vectors */
 44:   VecDuplicateVecs(x,3,&neP->workv);
 45:   w    = neP->workv[0];

 47:   /* Set components of vector w to random numbers */
 48:   PetscRandomCreate(PetscObjectComm((PetscObject)snes),&rctx);
 49:   PetscRandomSetFromOptions(rctx);
 50:   VecSetRandom(w,rctx);
 51:   PetscRandomDestroy(&rctx);

 53:   /* Open output file */
 54:   PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN,&flg);
 55:   if (flg) neP->fp = fopen(noise_file,"w");
 56:   else     neP->fp = fopen("noise.out","w");
 57:   if (!neP->fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file");
 58:   PetscInfo(snes,"Creating Jorge's differencing parameter context\n");

 60:   *outneP = neP;
 61:   return(0);
 62: }

 66: PetscErrorCode SNESDiffParameterDestroy_More(void *nePv)
 67: {
 68:   DIFFPAR_MORE   *neP = (DIFFPAR_MORE*)nePv;
 70:   int            err;

 73:   /* Destroy work vectors and close output file */
 74:   VecDestroyVecs(3,&neP->workv);
 75:   err  = fclose(neP->fp);
 76:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
 77:   PetscFree(neP);
 78:   return(0);
 79: }

 83: PetscErrorCode SNESDiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt)
 84: {
 85:   DIFFPAR_MORE   *neP = (DIFFPAR_MORE*)nePv;
 86:   Vec            w, xp, fvec;    /* work vectors to use in computing h */
 87:   double         zero = 0.0, hl, hu, h, fnoise_s, fder2_s;
 88:   PetscScalar    alpha;
 89:   PetscScalar    fval[7], tab[7][7], eps[7], f = -1;
 90:   double         rerrf = -1., fder2;
 92:   PetscInt       iter, k, i, j,  info;
 93:   PetscInt       nf = 7;         /* number of function evaluations */
 94:   PetscInt       fcount;
 95:   MPI_Comm       comm;
 96:   FILE           *fp;
 97:   PetscBool      noise_test = PETSC_FALSE;

100:   PetscObjectGetComm((PetscObject)snes,&comm);
101:   /* Call to SNESSetUp() just to set data structures in SNES context */
102:   if (!snes->setupcalled) {SNESSetUp(snes);}

104:   w    = neP->workv[0];
105:   xp   = neP->workv[1];
106:   fvec = neP->workv[2];
107:   fp   = neP->fp;

109:   /* Initialize parameters */
110:   hl       = zero;
111:   hu       = zero;
112:   h        = neP->h_first_try;
113:   fnoise_s = zero;
114:   fder2_s  = zero;
115:   fcount   = neP->function_count;

117:   /* We have 5 tries to attempt to compute a good hopt value */
118:   SNESGetIterationNumber(snes,&i);
119:   PetscFPrintf(comm,fp,"\n ------- SNES iteration %D ---------\n",i);
120:   for (iter=0; iter<5; iter++) {
121:     neP->h_first_try = h;

123:     /* Compute the nf function values needed to estimate the noise from
124:        the difference table */
125:     for (k=0; k<nf; k++) {
126:       alpha = h * (k+1 - (nf+1)/2);
127:       VecWAXPY(xp,alpha,p,x);
128:       SNESComputeFunction(snes,xp,fvec);
129:       neP->function_count++;
130:       VecDot(fvec,w,&fval[k]);
131:     }
132:     f = fval[(nf+1)/2 - 1];

134:     /* Construct the difference table */
135:     for (i=0; i<nf; i++) tab[i][0] = fval[i];

137:     for (j=0; j<6; j++) {
138:       for (i=0; i<nf-j; i++) {
139:         tab[i][j+1] = tab[i+1][j] - tab[i][j];
140:       }
141:     }

143:     /* Print the difference table */
144:     PetscFPrintf(comm,fp,"Difference Table: iter = %D\n",iter);
145:     for (i=0; i<nf; i++) {
146:       for (j=0; j<nf-i; j++) {
147:         PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);
148:       }
149:       PetscFPrintf(comm,fp,"\n");
150:     }

152:     /* Call the noise estimator */
153:     SNESNoise_dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);

155:     /* Output statements */
156:     rerrf = *fnoise/PetscAbsScalar(f);
157:     if (info == 1) {PetscFPrintf(comm,fp,"%s\n","Noise detected");}
158:     if (info == 2) {PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");}
159:     if (info == 3) {PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");}
160:     if (info == 4) {PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");}
161:     PetscFPrintf(comm,fp,"Approximate epsfcn %g  %g  %g  %g  %g  %g\n",(double)eps[0],(double)eps[1],(double)eps[2],(double)eps[3],(double)eps[4],(double)eps[5]);
162:     PetscFPrintf(comm,fp,"h = %g, fnoise = %g, fder2 = %g, rerrf = %g, hopt = %g\n\n",(double)h, (double)*fnoise, (double)fder2, (double)rerrf, (double)*hopt);

164:     /* Save fnoise and fder2. */
165:     if (*fnoise) fnoise_s = *fnoise;
166:     if (fder2) fder2_s = fder2;

168:     /* Check for noise detection. */
169:     if (fnoise_s && fder2_s) {
170:       *fnoise = fnoise_s;
171:       fder2   = fder2_s;
172:       *hopt   = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2));
173:       goto theend;
174:     } else {

176:       /* Update hl and hu, and determine new h */
177:       if (info == 2 || info == 4) {
178:         hl = h;
179:         if (hu == zero) h = 100*h;
180:         else            h = PetscMin(100*h,0.1*hu);
181:       } else if (info == 3) {
182:         hu = h;
183:         h  = PetscMax(1.0e-3,sqrt(hl/hu))*hu;
184:       }
185:     }
186:   }
187: theend:

189:   if (*fnoise < neP->fnoise_min) {
190:     PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",(double)*fnoise,(double)neP->fnoise_min);
191:     *fnoise = neP->fnoise_min;
192:     neP->fnoise_resets++;
193:   }
194:   if (*hopt < neP->hopt_min) {
195:     PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",(double)*hopt,(double)neP->hopt_min);
196:     *hopt = neP->hopt_min;
197:     neP->hopt_resets++;
198:   }

200:   PetscFPrintf(comm,fp,"Errors in derivative:\n");
201:   PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",(double)f,(double)*fnoise,(double)fder2,(double)*hopt);

203:   /* For now, compute h **each** MV Mult!! */
204:   /*
205:   PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg);
206:   if (!flg) {
207:     Mat mat;
208:     SNESGetJacobian(snes,&mat,NULL,NULL);
209:     SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);
210:   }
211:   */
212:   fcount = neP->function_count - fcount;
213:   PetscInfo5(snes,"fct_now = %D, fct_cum = %D, rerrf=%g, sqrt(noise)=%g, h_more=%g\n",fcount,neP->function_count,(double)rerrf,(double)PetscSqrtReal(*fnoise),(double)*hopt);

215:   PetscOptionsGetBool(NULL,"-noise_test",&noise_test,NULL);
216:   if (noise_test) {
217:     JacMatMultCompare(snes,x,p,*hopt);
218:   }
219:   return(0);
220: }

224: PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt)
225: {
226:   Vec            yy1, yy2; /* work vectors */
227:   PetscViewer    view2;    /* viewer */
228:   Mat            J;        /* analytic Jacobian (set as preconditioner matrix) */
229:   Mat            Jmf;      /* matrix-free Jacobian (set as true system matrix) */
230:   double         h;        /* differencing parameter */
231:   Vec            f;
232:   PetscScalar    alpha;
233:   PetscReal      yy1n,yy2n,enorm;
235:   PetscInt       i;
236:   PetscBool      printv = PETSC_FALSE;
237:   char           filename[32];
238:   MPI_Comm       comm;

241:   PetscObjectGetComm((PetscObject)snes,&comm);
242:   /* Compute function and analytic Jacobian at x */
243:   SNESGetJacobian(snes,&Jmf,&J,NULL,NULL);
244:   SNESComputeJacobian(snes,x,Jmf,J);
245:   SNESGetFunction(snes,&f,NULL,NULL);
246:   SNESComputeFunction(snes,x,f);

248:   /* Duplicate work vectors */
249:   VecDuplicate(x,&yy2);
250:   VecDuplicate(x,&yy1);

252:   /* Compute true matrix-vector product */
253:   MatMult(J,p,yy1);
254:   VecNorm(yy1,NORM_2,&yy1n);

256:   /* View product vector if desired */
257:   PetscOptionsGetBool(NULL,"-print_vecs",&printv,NULL);
258:   if (printv) {
259:     PetscViewerASCIIOpen(comm,"y1.out",&view2);
260:     PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);
261:     VecView(yy1,view2);
262:     PetscViewerDestroy(&view2);
263:   }

265:   /* Test Jacobian-vector product computation */
266:   alpha = -1.0;
267:   h     = 0.01 * hopt;
268:   for (i=0; i<5; i++) {
269:     /* Set differencing parameter for matrix-free multiplication */
270:     SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);

272:     /* Compute matrix-vector product via differencing approximation */
273:     MatMult(Jmf,p,yy2);
274:     VecNorm(yy2,NORM_2,&yy2n);

276:     /* View product vector if desired */
277:     if (printv) {
278:       sprintf(filename,"y2.%d.out",(int)i);
279:       PetscViewerASCIIOpen(comm,filename,&view2);
280:       PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);
281:       VecView(yy2,view2);
282:       PetscViewerDestroy(&view2);
283:     }

285:     /* Compute relative error */
286:     VecAXPY(yy2,alpha,yy1);
287:     VecNorm(yy2,NORM_2,&enorm);
288:     enorm = enorm/yy1n;
289:     PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",(double)h,(double)enorm);
290:     h    *= 10.0;
291:   }
292:   return(0);
293: }

295: static PetscInt lin_its_total = 0;

297: PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy)
298: {
300:   PetscInt       lin_its;

303:   SNESGetLinearSolveIterations(snes,&lin_its);
304:   lin_its_total += lin_its;
305:   PetscPrintf(PetscObjectComm((PetscObject)snes), "iter = %D, SNES Function norm = %g, lin_its = %D, total_lin_its = %D\n",its,(double)fnorm,lin_its,lin_its_total);

307:   SNESUnSetMatrixFreeParameter(snes);
308:   return(0);
309: }