13d177a5cSEmil Constantinescu 23d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I "petscts.h" I*/ 33d177a5cSEmil Constantinescu 43d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAdaptList; 53d177a5cSEmil Constantinescu static PetscBool TSGLLEAdaptPackageInitialized; 63d177a5cSEmil Constantinescu static PetscBool TSGLLEAdaptRegisterAllCalled; 73d177a5cSEmil Constantinescu static PetscClassId TSGLLEADAPT_CLASSID; 83d177a5cSEmil Constantinescu 93d177a5cSEmil Constantinescu struct _TSGLLEAdaptOps { 103d177a5cSEmil Constantinescu PetscErrorCode (*choose)(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool*); 113d177a5cSEmil Constantinescu PetscErrorCode (*destroy)(TSGLLEAdapt); 123d177a5cSEmil Constantinescu PetscErrorCode (*view)(TSGLLEAdapt,PetscViewer); 133d177a5cSEmil Constantinescu PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSGLLEAdapt); 143d177a5cSEmil Constantinescu }; 153d177a5cSEmil Constantinescu 163d177a5cSEmil Constantinescu struct _p_TSGLLEAdapt { 173d177a5cSEmil Constantinescu PETSCHEADER(struct _TSGLLEAdaptOps); 183d177a5cSEmil Constantinescu void *data; 193d177a5cSEmil Constantinescu }; 203d177a5cSEmil Constantinescu 213d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt); 223d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt); 233d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt); 243d177a5cSEmil Constantinescu 253d177a5cSEmil Constantinescu /*@C 263d177a5cSEmil Constantinescu TSGLLEAdaptRegister - adds a TSGLLEAdapt implementation 273d177a5cSEmil Constantinescu 283d177a5cSEmil Constantinescu Not Collective 293d177a5cSEmil Constantinescu 303d177a5cSEmil Constantinescu Input Parameters: 313d177a5cSEmil Constantinescu + name_scheme - name of user-defined adaptivity scheme 323d177a5cSEmil Constantinescu - routine_create - routine to create method context 333d177a5cSEmil Constantinescu 343d177a5cSEmil Constantinescu Notes: 353d177a5cSEmil Constantinescu TSGLLEAdaptRegister() may be called multiple times to add several user-defined families. 363d177a5cSEmil Constantinescu 373d177a5cSEmil Constantinescu Sample usage: 383d177a5cSEmil Constantinescu .vb 393d177a5cSEmil Constantinescu TSGLLEAdaptRegister("my_scheme",MySchemeCreate); 403d177a5cSEmil Constantinescu .ve 413d177a5cSEmil Constantinescu 423d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 433d177a5cSEmil Constantinescu $ TSGLLEAdaptSetType(ts,"my_scheme") 443d177a5cSEmil Constantinescu or at runtime via the option 453d177a5cSEmil Constantinescu $ -ts_adapt_type my_scheme 463d177a5cSEmil Constantinescu 473d177a5cSEmil Constantinescu Level: advanced 483d177a5cSEmil Constantinescu 493d177a5cSEmil Constantinescu .seealso: TSGLLEAdaptRegisterAll() 503d177a5cSEmil Constantinescu @*/ 513d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptRegister(const char sname[],PetscErrorCode (*function)(TSGLLEAdapt)) 523d177a5cSEmil Constantinescu { 533d177a5cSEmil Constantinescu PetscErrorCode ierr; 543d177a5cSEmil Constantinescu 553d177a5cSEmil Constantinescu PetscFunctionBegin; 561d36bdfdSBarry Smith ierr = TSGLLEAdaptInitializePackage();CHKERRQ(ierr); 573d177a5cSEmil Constantinescu ierr = PetscFunctionListAdd(&TSGLLEAdaptList,sname,function);CHKERRQ(ierr); 583d177a5cSEmil Constantinescu PetscFunctionReturn(0); 593d177a5cSEmil Constantinescu } 603d177a5cSEmil Constantinescu 613d177a5cSEmil Constantinescu /*@C 623d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in TSGLLEAdapt 633d177a5cSEmil Constantinescu 643d177a5cSEmil Constantinescu Not Collective 653d177a5cSEmil Constantinescu 663d177a5cSEmil Constantinescu Level: advanced 673d177a5cSEmil Constantinescu 683d177a5cSEmil Constantinescu .seealso: TSGLLEAdaptRegisterDestroy() 693d177a5cSEmil Constantinescu @*/ 703d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptRegisterAll(void) 713d177a5cSEmil Constantinescu { 723d177a5cSEmil Constantinescu PetscErrorCode ierr; 733d177a5cSEmil Constantinescu 743d177a5cSEmil Constantinescu PetscFunctionBegin; 753d177a5cSEmil Constantinescu if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(0); 763d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_TRUE; 773d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegister(TSGLLEADAPT_NONE,TSGLLEAdaptCreate_None);CHKERRQ(ierr); 783d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegister(TSGLLEADAPT_SIZE,TSGLLEAdaptCreate_Size);CHKERRQ(ierr); 793d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegister(TSGLLEADAPT_BOTH,TSGLLEAdaptCreate_Both);CHKERRQ(ierr); 803d177a5cSEmil Constantinescu PetscFunctionReturn(0); 813d177a5cSEmil Constantinescu } 823d177a5cSEmil Constantinescu 833d177a5cSEmil Constantinescu /*@C 843d177a5cSEmil Constantinescu TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 853d177a5cSEmil Constantinescu called from PetscFinalize(). 863d177a5cSEmil Constantinescu 873d177a5cSEmil Constantinescu Level: developer 883d177a5cSEmil Constantinescu 893d177a5cSEmil Constantinescu .seealso: PetscFinalize() 903d177a5cSEmil Constantinescu @*/ 913d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptFinalizePackage(void) 923d177a5cSEmil Constantinescu { 933d177a5cSEmil Constantinescu PetscErrorCode ierr; 943d177a5cSEmil Constantinescu 953d177a5cSEmil Constantinescu PetscFunctionBegin; 963d177a5cSEmil Constantinescu ierr = PetscFunctionListDestroy(&TSGLLEAdaptList);CHKERRQ(ierr); 973d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_FALSE; 983d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_FALSE; 993d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1003d177a5cSEmil Constantinescu } 1013d177a5cSEmil Constantinescu 1023d177a5cSEmil Constantinescu /*@C 1033d177a5cSEmil Constantinescu TSGLLEAdaptInitializePackage - This function initializes everything in the TSGLLEAdapt package. It is 1048a690491SBarry Smith called from TSInitializePackage(). 1053d177a5cSEmil Constantinescu 1063d177a5cSEmil Constantinescu Level: developer 1073d177a5cSEmil Constantinescu 1083d177a5cSEmil Constantinescu .seealso: PetscInitialize() 1093d177a5cSEmil Constantinescu @*/ 1103d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptInitializePackage(void) 1113d177a5cSEmil Constantinescu { 1123d177a5cSEmil Constantinescu PetscErrorCode ierr; 1133d177a5cSEmil Constantinescu 1143d177a5cSEmil Constantinescu PetscFunctionBegin; 1153d177a5cSEmil Constantinescu if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(0); 1163d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_TRUE; 1173d177a5cSEmil Constantinescu ierr = PetscClassIdRegister("TSGLLEAdapt",&TSGLLEADAPT_CLASSID);CHKERRQ(ierr); 1183d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegisterAll();CHKERRQ(ierr); 1193d177a5cSEmil Constantinescu ierr = PetscRegisterFinalize(TSGLLEAdaptFinalizePackage);CHKERRQ(ierr); 1203d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1213d177a5cSEmil Constantinescu } 1223d177a5cSEmil Constantinescu 1233d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt,TSGLLEAdaptType type) 1243d177a5cSEmil Constantinescu { 1253d177a5cSEmil Constantinescu PetscErrorCode ierr,(*r)(TSGLLEAdapt); 1263d177a5cSEmil Constantinescu 1273d177a5cSEmil Constantinescu PetscFunctionBegin; 1283d177a5cSEmil Constantinescu ierr = PetscFunctionListFind(TSGLLEAdaptList,type,&r);CHKERRQ(ierr); 1293d177a5cSEmil Constantinescu if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAdapt type \"%s\" given",type); 1303d177a5cSEmil Constantinescu if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1313d177a5cSEmil Constantinescu ierr = (*r)(adapt);CHKERRQ(ierr); 1323d177a5cSEmil Constantinescu ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 1333d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1343d177a5cSEmil Constantinescu } 1353d177a5cSEmil Constantinescu 1363d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt,const char prefix[]) 1373d177a5cSEmil Constantinescu { 1383d177a5cSEmil Constantinescu PetscErrorCode ierr; 1393d177a5cSEmil Constantinescu 1403d177a5cSEmil Constantinescu PetscFunctionBegin; 1413d177a5cSEmil Constantinescu ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 1423d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1433d177a5cSEmil Constantinescu } 1443d177a5cSEmil Constantinescu 1453d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt,PetscViewer viewer) 1463d177a5cSEmil Constantinescu { 1473d177a5cSEmil Constantinescu PetscErrorCode ierr; 1483d177a5cSEmil Constantinescu PetscBool iascii; 1493d177a5cSEmil Constantinescu 1503d177a5cSEmil Constantinescu PetscFunctionBegin; 1513d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1523d177a5cSEmil Constantinescu if (iascii) { 1533d177a5cSEmil Constantinescu ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 1543d177a5cSEmil Constantinescu if (adapt->ops->view) { 1553d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1563d177a5cSEmil Constantinescu ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 1573d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1583d177a5cSEmil Constantinescu } 1593d177a5cSEmil Constantinescu } 1603d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1613d177a5cSEmil Constantinescu } 1623d177a5cSEmil Constantinescu 1633d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt) 1643d177a5cSEmil Constantinescu { 1653d177a5cSEmil Constantinescu PetscErrorCode ierr; 1663d177a5cSEmil Constantinescu 1673d177a5cSEmil Constantinescu PetscFunctionBegin; 1683d177a5cSEmil Constantinescu if (!*adapt) PetscFunctionReturn(0); 1693d177a5cSEmil Constantinescu PetscValidHeaderSpecific(*adapt,TSGLLEADAPT_CLASSID,1); 170*c793f718SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 1713d177a5cSEmil Constantinescu if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 1723d177a5cSEmil Constantinescu ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 1733d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1743d177a5cSEmil Constantinescu } 1753d177a5cSEmil Constantinescu 1763d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSGLLEAdapt adapt) 1773d177a5cSEmil Constantinescu { 1783d177a5cSEmil Constantinescu PetscErrorCode ierr; 1793d177a5cSEmil Constantinescu char type[256] = TSGLLEADAPT_BOTH; 1803d177a5cSEmil Constantinescu PetscBool flg; 1813d177a5cSEmil Constantinescu 1823d177a5cSEmil Constantinescu PetscFunctionBegin; 1833d177a5cSEmil Constantinescu /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this 1843d177a5cSEmil Constantinescu * function can only be called from inside TSSetFromOptions_GLLE() */ 1853d177a5cSEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"TSGLLE Adaptivity options");CHKERRQ(ierr); 1863d177a5cSEmil Constantinescu ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSGLLEAdaptSetType",TSGLLEAdaptList, 1873d177a5cSEmil Constantinescu ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 1883d177a5cSEmil Constantinescu if (flg || !((PetscObject)adapt)->type_name) { 1893d177a5cSEmil Constantinescu ierr = TSGLLEAdaptSetType(adapt,type);CHKERRQ(ierr); 1903d177a5cSEmil Constantinescu } 1913d177a5cSEmil Constantinescu if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 1923d177a5cSEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 1933d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1943d177a5cSEmil Constantinescu } 1953d177a5cSEmil Constantinescu 1963d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool *finish) 1973d177a5cSEmil Constantinescu { 1983d177a5cSEmil Constantinescu PetscErrorCode ierr; 1993d177a5cSEmil Constantinescu 2003d177a5cSEmil Constantinescu PetscFunctionBegin; 2013d177a5cSEmil Constantinescu PetscValidHeaderSpecific(adapt,TSGLLEADAPT_CLASSID,1); 2023d177a5cSEmil Constantinescu PetscValidIntPointer(orders,3); 2033d177a5cSEmil Constantinescu PetscValidPointer(errors,4); 2043d177a5cSEmil Constantinescu PetscValidPointer(cost,5); 2053d177a5cSEmil Constantinescu PetscValidIntPointer(next_sc,9); 2063d177a5cSEmil Constantinescu PetscValidPointer(next_h,10); 2073d177a5cSEmil Constantinescu PetscValidIntPointer(finish,11); 2083d177a5cSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,n,orders,errors,cost,cur,h,tleft,next_sc,next_h,finish);CHKERRQ(ierr); 2093d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2103d177a5cSEmil Constantinescu } 2113d177a5cSEmil Constantinescu 2123d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm,TSGLLEAdapt *inadapt) 2133d177a5cSEmil Constantinescu { 2143d177a5cSEmil Constantinescu PetscErrorCode ierr; 2153d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 2163d177a5cSEmil Constantinescu 2173d177a5cSEmil Constantinescu PetscFunctionBegin; 2183d177a5cSEmil Constantinescu *inadapt = NULL; 2193d177a5cSEmil Constantinescu ierr = PetscHeaderCreate(adapt,TSGLLEADAPT_CLASSID,"TSGLLEAdapt","General Linear adaptivity","TS",comm,TSGLLEAdaptDestroy,TSGLLEAdaptView);CHKERRQ(ierr); 2203d177a5cSEmil Constantinescu *inadapt = adapt; 2213d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2223d177a5cSEmil Constantinescu } 2233d177a5cSEmil Constantinescu 2243d177a5cSEmil Constantinescu 2253d177a5cSEmil Constantinescu /* 2263d177a5cSEmil Constantinescu * Implementations 2273d177a5cSEmil Constantinescu */ 2283d177a5cSEmil Constantinescu 2293d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt) 2303d177a5cSEmil Constantinescu { 2313d177a5cSEmil Constantinescu PetscErrorCode ierr; 2323d177a5cSEmil Constantinescu 2333d177a5cSEmil Constantinescu PetscFunctionBegin; 2343d177a5cSEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 2353d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2363d177a5cSEmil Constantinescu } 2373d177a5cSEmil Constantinescu 2383d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */ 2393d177a5cSEmil Constantinescu typedef struct { 2403d177a5cSEmil Constantinescu PetscInt scheme; 2413d177a5cSEmil Constantinescu PetscReal h; 2423d177a5cSEmil Constantinescu } TSGLLEAdapt_None; 2433d177a5cSEmil Constantinescu 2443d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAdaptChoose_None(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool *finish) 2453d177a5cSEmil Constantinescu { 2463d177a5cSEmil Constantinescu 2473d177a5cSEmil Constantinescu PetscFunctionBegin; 2483d177a5cSEmil Constantinescu *next_sc = cur; 2493d177a5cSEmil Constantinescu *next_h = h; 2503d177a5cSEmil Constantinescu if (*next_h > tleft) { 2513d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2523d177a5cSEmil Constantinescu *next_h = tleft; 2533d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2543d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2553d177a5cSEmil Constantinescu } 2563d177a5cSEmil Constantinescu 2573d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt) 2583d177a5cSEmil Constantinescu { 2593d177a5cSEmil Constantinescu PetscErrorCode ierr; 2603d177a5cSEmil Constantinescu TSGLLEAdapt_None *a; 2613d177a5cSEmil Constantinescu 2623d177a5cSEmil Constantinescu PetscFunctionBegin; 2633d177a5cSEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 2643d177a5cSEmil Constantinescu adapt->data = (void*)a; 2653d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_None; 2663d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2673d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2683d177a5cSEmil Constantinescu } 2693d177a5cSEmil Constantinescu 2703d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */ 2713d177a5cSEmil Constantinescu typedef struct { 2723d177a5cSEmil Constantinescu PetscReal desired_h; 2733d177a5cSEmil Constantinescu } TSGLLEAdapt_Size; 2743d177a5cSEmil Constantinescu 2753d177a5cSEmil Constantinescu 2763d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAdaptChoose_Size(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool *finish) 2773d177a5cSEmil Constantinescu { 2783d177a5cSEmil Constantinescu TSGLLEAdapt_Size *sz = (TSGLLEAdapt_Size*)adapt->data; 2793d177a5cSEmil Constantinescu PetscReal dec = 0.2,inc = 5.0,safe = 0.9,optimal,last_desired_h; 2803d177a5cSEmil Constantinescu 2813d177a5cSEmil Constantinescu PetscFunctionBegin; 2823d177a5cSEmil Constantinescu *next_sc = cur; 2833d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[cur],(PetscReal)-1./(safe*orders[cur])); 2843d177a5cSEmil Constantinescu /* Step sizes oscillate when there is no smoothing. Here we use a geometric mean of the current step size and the 2853d177a5cSEmil Constantinescu * one that would have been taken (without smoothing) on the last step. */ 2863d177a5cSEmil Constantinescu last_desired_h = sz->desired_h; 2873d177a5cSEmil Constantinescu sz->desired_h = h*PetscMax(dec,PetscMin(inc,optimal)); /* Trim to [dec,inc] */ 2883d177a5cSEmil Constantinescu 2893d177a5cSEmil Constantinescu /* Normally only happens on the first step */ 2903d177a5cSEmil Constantinescu if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h); 2913d177a5cSEmil Constantinescu else *next_h = sz->desired_h; 2923d177a5cSEmil Constantinescu 2933d177a5cSEmil Constantinescu if (*next_h > tleft) { 2943d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2953d177a5cSEmil Constantinescu *next_h = tleft; 2963d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2973d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2983d177a5cSEmil Constantinescu } 2993d177a5cSEmil Constantinescu 3003d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt) 3013d177a5cSEmil Constantinescu { 3023d177a5cSEmil Constantinescu PetscErrorCode ierr; 3033d177a5cSEmil Constantinescu TSGLLEAdapt_Size *a; 3043d177a5cSEmil Constantinescu 3053d177a5cSEmil Constantinescu PetscFunctionBegin; 3063d177a5cSEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 3073d177a5cSEmil Constantinescu adapt->data = (void*)a; 3083d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Size; 3093d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 3103d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3113d177a5cSEmil Constantinescu } 3123d177a5cSEmil Constantinescu 3133d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */ 3143d177a5cSEmil Constantinescu typedef struct { 3153d177a5cSEmil Constantinescu PetscInt count_at_order; 3163d177a5cSEmil Constantinescu PetscReal desired_h; 3173d177a5cSEmil Constantinescu } TSGLLEAdapt_Both; 3183d177a5cSEmil Constantinescu 3193d177a5cSEmil Constantinescu 3203d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAdaptChoose_Both(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool *finish) 3213d177a5cSEmil Constantinescu { 3223d177a5cSEmil Constantinescu TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both*)adapt->data; 3233d177a5cSEmil Constantinescu PetscErrorCode ierr; 3243d177a5cSEmil Constantinescu PetscReal dec = 0.2,inc = 5.0,safe = 0.9; 3253d177a5cSEmil Constantinescu struct {PetscInt id; PetscReal h,eff;} best={-1,0,0},trial={-1,0,0},current={-1,0,0}; 3263d177a5cSEmil Constantinescu PetscInt i; 3273d177a5cSEmil Constantinescu 3283d177a5cSEmil Constantinescu PetscFunctionBegin; 3293d177a5cSEmil Constantinescu for (i=0; i<n; i++) { 3303d177a5cSEmil Constantinescu PetscReal optimal; 3313d177a5cSEmil Constantinescu trial.id = i; 3323d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[i],(PetscReal)-1./(safe*orders[i])); 3333d177a5cSEmil Constantinescu trial.h = h*optimal; 3343d177a5cSEmil Constantinescu trial.eff = trial.h/cost[i]; 335580bdb30SBarry Smith if (trial.eff > best.eff) {ierr = PetscArraycpy(&best,&trial,1);CHKERRQ(ierr);} 336580bdb30SBarry Smith if (i == cur) {ierr = PetscArraycpy(¤t,&trial,1);CHKERRQ(ierr);} 3373d177a5cSEmil Constantinescu } 3383d177a5cSEmil Constantinescu /* Only switch orders if the scheme offers significant benefits over the current one. 3393d177a5cSEmil Constantinescu When the scheme is not changing, only change step size if it offers significant benefits. */ 3403d177a5cSEmil Constantinescu if (best.eff < 1.2*current.eff || both->count_at_order < orders[cur]+2) { 3413d177a5cSEmil Constantinescu PetscReal last_desired_h; 3423d177a5cSEmil Constantinescu *next_sc = current.id; 3433d177a5cSEmil Constantinescu last_desired_h = both->desired_h; 3443d177a5cSEmil Constantinescu both->desired_h = PetscMax(h*dec,PetscMin(h*inc,current.h)); 3453d177a5cSEmil Constantinescu *next_h = (both->count_at_order > 0) 3463d177a5cSEmil Constantinescu ? PetscSqrtReal(last_desired_h * both->desired_h) 3473d177a5cSEmil Constantinescu : both->desired_h; 3483d177a5cSEmil Constantinescu both->count_at_order++; 3493d177a5cSEmil Constantinescu } else { 3503d177a5cSEmil Constantinescu PetscReal rat = cost[best.id]/cost[cur]; 3513d177a5cSEmil Constantinescu *next_sc = best.id; 3523d177a5cSEmil Constantinescu *next_h = PetscMax(h*rat*dec,PetscMin(h*rat*inc,best.h)); 3533d177a5cSEmil Constantinescu both->count_at_order = 0; 3543d177a5cSEmil Constantinescu both->desired_h = best.h; 3553d177a5cSEmil Constantinescu } 3563d177a5cSEmil Constantinescu 3573d177a5cSEmil Constantinescu if (*next_h > tleft) { 3583d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 3593d177a5cSEmil Constantinescu *next_h = tleft; 3603d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 3613d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3623d177a5cSEmil Constantinescu } 3633d177a5cSEmil Constantinescu 3643d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt) 3653d177a5cSEmil Constantinescu { 3663d177a5cSEmil Constantinescu PetscErrorCode ierr; 3673d177a5cSEmil Constantinescu TSGLLEAdapt_Both *a; 3683d177a5cSEmil Constantinescu 3693d177a5cSEmil Constantinescu PetscFunctionBegin; 3703d177a5cSEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 3713d177a5cSEmil Constantinescu adapt->data = (void*)a; 3723d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Both; 3733d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 3743d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3753d177a5cSEmil Constantinescu } 376