Actual source code: gladapt.c

  1: #define PETSCTS_DLL

 3:  #include gl.h

  5: static PetscFList TSGLAdaptList;
  6: static PetscTruth TSGLAdaptPackageInitialized;
  7: static PetscTruth TSGLAdaptRegisterAllCalled;
  8: static PetscCookie TSGLADAPT_COOKIE;

 10: struct _TSGLAdaptOps {
 11:   PetscErrorCode (*choose)(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscTruth*);
 12:   PetscErrorCode (*destroy)(TSGLAdapt);
 13:   PetscErrorCode (*view)(TSGLAdapt,PetscViewer);
 14:   PetscErrorCode (*setfromoptions)(TSGLAdapt);
 15: };

 17: struct _p_TSGLAdapt {
 18:   PETSCHEADER(struct _TSGLAdaptOps);
 19:   void *data;
 20: };

 22: static PetscErrorCode TSGLAdaptCreate_None(TSGLAdapt);
 23: static PetscErrorCode TSGLAdaptCreate_Size(TSGLAdapt);
 24: static PetscErrorCode TSGLAdaptCreate_Both(TSGLAdapt);

 28: /*@C
 29:    TSGLAdaptRegister - see TSGLAdaptRegisterDynamic()

 31:    Level: advanced
 32: @*/
 33: PetscErrorCode  TSGLAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSGLAdapt))
 34: {
 36:   char           fullname[PETSC_MAX_PATH_LEN];

 39:   PetscFListConcat(path,name,fullname);
 40:   PetscFListAdd(&TSGLAdaptList,sname,fullname,(void(*)(void))function);
 41:   return(0);
 42: }

 46: /*@C
 47:   TSGLAdaptRegisterAll - Registers all of the adaptivity schemes in TSGLAdapt

 49:   Not Collective

 51:   Level: advanced

 53: .keywords: TSGLAdapt, register, all

 55: .seealso: TSGLAdaptRegisterDestroy()
 56: @*/
 57: PetscErrorCode  TSGLAdaptRegisterAll(const char path[])
 58: {

 62:   TSGLAdaptRegisterDynamic(TSGLADAPT_NONE,path,"TSGLAdaptCreate_None",TSGLAdaptCreate_None);
 63:   TSGLAdaptRegisterDynamic(TSGLADAPT_SIZE,path,"TSGLAdaptCreate_Size",TSGLAdaptCreate_Size);
 64:   TSGLAdaptRegisterDynamic(TSGLADAPT_BOTH,path,"TSGLAdaptCreate_Both",TSGLAdaptCreate_Both);
 65:   return(0);
 66: }

 70: /*@C
 71:   TSGLFinalizePackage - This function destroys everything in the TSGL package. It is
 72:   called from PetscFinalize().

 74:   Level: developer

 76: .keywords: Petsc, destroy, package
 77: .seealso: PetscFinalize()
 78: @*/
 79: PetscErrorCode  TSGLAdaptFinalizePackage(void)
 80: {
 82:   TSGLAdaptPackageInitialized = PETSC_FALSE;
 83:   TSGLAdaptRegisterAllCalled  = PETSC_FALSE;
 84:   TSGLAdaptList               = PETSC_NULL;
 85:   return(0);
 86: }

 90: /*@C
 91:   TSGLAdaptInitializePackage - This function initializes everything in the TSGLAdapt package. It is
 92:   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
 93:   TSCreate_GL() when using static libraries.

 95:   Input Parameter:
 96:   path - The dynamic library path, or PETSC_NULL

 98:   Level: developer

100: .keywords: TSGLAdapt, initialize, package
101: .seealso: PetscInitialize()
102: @*/
103: PetscErrorCode  TSGLAdaptInitializePackage(const char path[])
104: {

108:   if (TSGLAdaptPackageInitialized) return(0);
109:   TSGLAdaptPackageInitialized = PETSC_TRUE;
110:   PetscCookieRegister("TSGLAdapt",&TSGLADAPT_COOKIE);
111:   TSGLAdaptRegisterAll(path);
112:   PetscRegisterFinalize(TSGLAdaptFinalizePackage);
113:   return(0);
114: }

118: /*@C
119:    TSGLAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSGLAdaptRegister()/TSGLAdaptRegisterDynamic().

121:    Not Collective

123:    Level: advanced

125: .keywords: TSGLAdapt, register, destroy
126: .seealso: TSGLAdaptRegister(), TSGLAdaptRegisterAll(), TSGLAdaptRegisterDynamic()
127: @*/
128: PetscErrorCode  TSGLAdaptRegisterDestroy(void)
129: {

133:   PetscFListDestroy(&TSGLAdaptList);
134:   TSGLAdaptRegisterAllCalled = PETSC_FALSE;
135:   return(0);
136: }


141: PetscErrorCode  TSGLAdaptSetType(TSGLAdapt adapt,const TSGLAdaptType type)
142: {
143:   PetscErrorCode ierr,(*r)(TSGLAdapt);

146:   PetscFListFind(TSGLAdaptList,((PetscObject)adapt)->comm,type,(void(**)(void))&r);
147:   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLAdapt type \"%s\" given",type);
148:   if (((PetscObject)adapt)->type_name) {(*adapt->ops->destroy)(adapt);}
149:   (*r)(adapt);
150:   PetscObjectChangeTypeName((PetscObject)adapt,type);
151:   return(0);
152: }

156: PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt adapt,const char prefix[])
157: {

161:   PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);
162:   return(0);
163: }

167: PetscErrorCode  TSGLAdaptView(TSGLAdapt adapt,PetscViewer viewer)
168: {
170:   PetscTruth iascii;

173:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
174:   if (iascii) {
175:     if (((PetscObject)adapt)->prefix) {
176:       PetscViewerASCIIPrintf(viewer,"TSGLAdapt object: (%s)\n",((PetscObject)adapt)->prefix);
177:     } else {
178:       PetscViewerASCIIPrintf(viewer,"TSGLAdapt object:\n");
179:     }
180:     PetscViewerASCIIPrintf(viewer,"  type: %s\n",((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:"(not set yet)");
181:     if (adapt->ops->view) {
182:       PetscViewerASCIIPushTab(viewer);
183:       (*adapt->ops->view)(adapt,viewer);
184:       PetscViewerASCIIPopTab(viewer);
185:     }
186:   } else {
187:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_GL",((PetscObject)viewer)->type_name);
188:   }

190:   return(0);
191: }

195: PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt adapt)
196: {

201:   if (--((PetscObject)adapt)->refct > 0) return(0);
202:   if (adapt->ops->destroy) {(*adapt->ops->destroy)(adapt);}
203:   PetscHeaderDestroy(adapt);
204:   return(0);
205: }

209: PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt adapt)
210: {
212:   char           type[256] = TSGLADAPT_BOTH;
213:   PetscTruth     flg;

216:   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGL, but currently this
217:   * function can only be called from inside TSSetFromOptions_GL()  */
218:   PetscOptionsHead("TSGL Adaptivity options");
219:   PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSGLAdaptSetType",TSGLAdaptList,
220:                           ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof type,&flg);
221:   if (flg || !((PetscObject)adapt)->type_name) {
222:     TSGLAdaptSetType(adapt,type);
223:   }
224:   if (adapt->ops->setfromoptions) {(*adapt->ops->setfromoptions)(adapt);}
225:   PetscOptionsTail();
226:   return(0);
227: }

231: PetscErrorCode  TSGLAdaptChoose(TSGLAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscTruth *finish)
232: {

243:   (*adapt->ops->choose)(adapt,n,orders,errors,cost,cur,h,tleft,next_sc,next_h,finish);
244:   return(0);
245: }

249: PetscErrorCode  TSGLAdaptCreate(MPI_Comm comm,TSGLAdapt *inadapt)
250: {
252:   TSGLAdapt adapt;

255:   *inadapt = 0;
256:   PetscHeaderCreate(adapt,_p_TSGLAdapt,struct _TSGLAdaptOps,TSGLADAPT_COOKIE,0,"TSGLAdapt",comm,TSGLAdaptDestroy,TSGLAdaptView);
257:   *inadapt = adapt;
258:   return(0);
259: }


262: /*
263: *  Implementations
264: */

268: static PetscErrorCode TSGLAdaptDestroy_JustFree(TSGLAdapt adapt)
269: {

273:   PetscFree(adapt->data);
274:   return(0);
275: }

277: /* -------------------------------- None ----------------------------------- */
278: typedef struct {
279:   PetscInt scheme;
280:   PetscReal h;
281: } TSGLAdapt_None;

285: static PetscErrorCode TSGLAdaptChoose_None(TSGLAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscTruth *finish)
286: {

289:   *next_sc = cur;
290:   *next_h = h;
291:   if (*next_h > tleft) {
292:     *finish = PETSC_TRUE;
293:     *next_h = tleft;
294:   } else {
295:     *finish = PETSC_FALSE;
296:   }
297:   return(0);
298: }

302: static PetscErrorCode TSGLAdaptCreate_None(TSGLAdapt adapt)
303: {
305:   TSGLAdapt_None *a;

308:   PetscNewLog(adapt,TSGLAdapt_None,&a);
309:   adapt->data = (void*)a;
310:   adapt->ops->choose = TSGLAdaptChoose_None;
311:   adapt->ops->destroy = TSGLAdaptDestroy_JustFree;
312:   return(0);
313: }


316: /* -------------------------------- Size ----------------------------------- */
317: typedef struct {
318:   PetscReal desired_h;
319: } TSGLAdapt_Size;


324: static PetscErrorCode TSGLAdaptChoose_Size(TSGLAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscTruth *finish)
325: {
326:   TSGLAdapt_Size *sz = (TSGLAdapt_Size*)adapt->data;
327:   PetscReal dec = 0.2,inc = 5.0,safe = 0.9,optimal,last_desired_h;

330:   *next_sc = cur;
331:   optimal = pow(errors[cur],-1./(safe*orders[cur]));
332:   /* Step sizes oscillate when there is no smoothing.  Here we use a geometric mean of the the current step size and the
333:   * one that would have been taken (without smoothing) on the last step. */
334:   last_desired_h = sz->desired_h;
335:   sz->desired_h = h*PetscMax(dec,PetscMin(inc,optimal)); /* Trim to [dec,inc] */
336:   if (last_desired_h > 1e-14) {                          /* Normally only happens on the first step */
337:     *next_h = sqrt(last_desired_h * sz->desired_h);
338:   } else {
339:     *next_h = sz->desired_h;
340:   }
341:   if (*next_h > tleft) {
342:     *finish = PETSC_TRUE;
343:     *next_h = tleft;
344:   } else {
345:     *finish = PETSC_FALSE;
346:   }
347:   return(0);
348: }

352: static PetscErrorCode TSGLAdaptCreate_Size(TSGLAdapt adapt)
353: {
355:   TSGLAdapt_Size *a;

358:   PetscNewLog(adapt,TSGLAdapt_Size,&a);
359:   adapt->data = (void*)a;
360:   adapt->ops->choose = TSGLAdaptChoose_Size;
361:   adapt->ops->destroy = TSGLAdaptDestroy_JustFree;
362:   return(0);
363: }

365: /* -------------------------------- Both ----------------------------------- */
366: typedef struct {
367:   PetscInt count_at_order;
368:   PetscReal desired_h;
369: } TSGLAdapt_Both;


374: static PetscErrorCode TSGLAdaptChoose_Both(TSGLAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscTruth *finish)
375: {
376:   TSGLAdapt_Both *both = (TSGLAdapt_Both*)adapt->data;
378:   PetscReal dec = 0.2,inc = 5.0,safe = 0.9;
379:   struct {PetscInt id; PetscReal h,eff;} best={-1,0,0},trial={-1,0,0},current={-1,0,0};
380:   PetscInt i;

383:   for (i=0; i<n; i++) {
384:     PetscReal optimal;
385:     trial.id = i;
386:     optimal = pow(errors[i],-1./(safe*orders[i]));
387:     trial.h = h*optimal;
388:     trial.eff = trial.h/cost[i];
389:     if (trial.eff > best.eff) {PetscMemcpy(&best,&trial,sizeof(trial));}
390:     if (i == cur) {PetscMemcpy(&current,&trial,sizeof(trial));}
391:   }
392:   /* Only switch orders if the scheme offers significant benefits over the current one.
393:   When the scheme is not changing, only change step size if it offers significant benefits. */
394:   if (best.eff < 1.2*current.eff || both->count_at_order < orders[cur]+2) {
395:     PetscReal last_desired_h;
396:     *next_sc = current.id;
397:     last_desired_h = both->desired_h;
398:     both->desired_h = PetscMax(h*dec,PetscMin(h*inc,current.h));
399:     *next_h = (both->count_at_order > 0)
400:       ? sqrt(last_desired_h * both->desired_h)
401:       : both->desired_h;
402:     both->count_at_order++;
403:   } else {
404:     PetscReal rat = cost[best.id]/cost[cur];
405:     *next_sc = best.id;
406:     *next_h  = PetscMax(h*rat*dec,PetscMin(h*rat*inc,best.h));
407:     both->count_at_order = 0;
408:     both->desired_h = best.h;
409:   }

411:   if (*next_h > tleft) {
412:     *finish = PETSC_TRUE;
413:     *next_h = tleft;
414:   } else {
415:     *finish = PETSC_FALSE;
416:   }
417:   return(0);
418: }

422: static PetscErrorCode TSGLAdaptCreate_Both(TSGLAdapt adapt)
423: {
425:   TSGLAdapt_Both *a;

428:   PetscNewLog(adapt,TSGLAdapt_Both,&a);
429:   adapt->data = (void*)a;
430:   adapt->ops->choose = TSGLAdaptChoose_Both;
431:   adapt->ops->destroy = TSGLAdaptDestroy_JustFree;
432:   return(0);
433: }