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 PetscFunctionBegin; 545f80ce2aSJacob Faibussowitsch CHKERRQ(TSGLLEAdaptInitializePackage()); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&TSGLLEAdaptList,sname,function)); 563d177a5cSEmil Constantinescu PetscFunctionReturn(0); 573d177a5cSEmil Constantinescu } 583d177a5cSEmil Constantinescu 593d177a5cSEmil Constantinescu /*@C 603d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in TSGLLEAdapt 613d177a5cSEmil Constantinescu 623d177a5cSEmil Constantinescu Not Collective 633d177a5cSEmil Constantinescu 643d177a5cSEmil Constantinescu Level: advanced 653d177a5cSEmil Constantinescu 663d177a5cSEmil Constantinescu .seealso: TSGLLEAdaptRegisterDestroy() 673d177a5cSEmil Constantinescu @*/ 683d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptRegisterAll(void) 693d177a5cSEmil Constantinescu { 703d177a5cSEmil Constantinescu PetscFunctionBegin; 713d177a5cSEmil Constantinescu if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(0); 723d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_TRUE; 735f80ce2aSJacob Faibussowitsch CHKERRQ(TSGLLEAdaptRegister(TSGLLEADAPT_NONE,TSGLLEAdaptCreate_None)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(TSGLLEAdaptRegister(TSGLLEADAPT_SIZE,TSGLLEAdaptCreate_Size)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(TSGLLEAdaptRegister(TSGLLEADAPT_BOTH,TSGLLEAdaptCreate_Both)); 763d177a5cSEmil Constantinescu PetscFunctionReturn(0); 773d177a5cSEmil Constantinescu } 783d177a5cSEmil Constantinescu 793d177a5cSEmil Constantinescu /*@C 803d177a5cSEmil Constantinescu TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 813d177a5cSEmil Constantinescu called from PetscFinalize(). 823d177a5cSEmil Constantinescu 833d177a5cSEmil Constantinescu Level: developer 843d177a5cSEmil Constantinescu 853d177a5cSEmil Constantinescu .seealso: PetscFinalize() 863d177a5cSEmil Constantinescu @*/ 873d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptFinalizePackage(void) 883d177a5cSEmil Constantinescu { 893d177a5cSEmil Constantinescu PetscFunctionBegin; 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListDestroy(&TSGLLEAdaptList)); 913d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_FALSE; 923d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_FALSE; 933d177a5cSEmil Constantinescu PetscFunctionReturn(0); 943d177a5cSEmil Constantinescu } 953d177a5cSEmil Constantinescu 963d177a5cSEmil Constantinescu /*@C 973d177a5cSEmil Constantinescu TSGLLEAdaptInitializePackage - This function initializes everything in the TSGLLEAdapt package. It is 988a690491SBarry Smith called from TSInitializePackage(). 993d177a5cSEmil Constantinescu 1003d177a5cSEmil Constantinescu Level: developer 1013d177a5cSEmil Constantinescu 1023d177a5cSEmil Constantinescu .seealso: PetscInitialize() 1033d177a5cSEmil Constantinescu @*/ 1043d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptInitializePackage(void) 1053d177a5cSEmil Constantinescu { 1063d177a5cSEmil Constantinescu PetscFunctionBegin; 1073d177a5cSEmil Constantinescu if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(0); 1083d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_TRUE; 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscClassIdRegister("TSGLLEAdapt",&TSGLLEADAPT_CLASSID)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(TSGLLEAdaptRegisterAll()); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRegisterFinalize(TSGLLEAdaptFinalizePackage)); 1123d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1133d177a5cSEmil Constantinescu } 1143d177a5cSEmil Constantinescu 1153d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt,TSGLLEAdaptType type) 1163d177a5cSEmil Constantinescu { 1175f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSGLLEAdapt); 1183d177a5cSEmil Constantinescu 1193d177a5cSEmil Constantinescu PetscFunctionBegin; 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListFind(TSGLLEAdaptList,type,&r)); 1213c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAdapt type \"%s\" given",type); 1225f80ce2aSJacob Faibussowitsch if (((PetscObject)adapt)->type_name) CHKERRQ((*adapt->ops->destroy)(adapt)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ((*r)(adapt)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)adapt,type)); 1253d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1263d177a5cSEmil Constantinescu } 1273d177a5cSEmil Constantinescu 1283d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt,const char prefix[]) 1293d177a5cSEmil Constantinescu { 1303d177a5cSEmil Constantinescu PetscFunctionBegin; 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix)); 1323d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1333d177a5cSEmil Constantinescu } 1343d177a5cSEmil Constantinescu 1353d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt,PetscViewer viewer) 1363d177a5cSEmil Constantinescu { 1373d177a5cSEmil Constantinescu PetscBool iascii; 1383d177a5cSEmil Constantinescu 1393d177a5cSEmil Constantinescu PetscFunctionBegin; 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 1413d177a5cSEmil Constantinescu if (iascii) { 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer)); 1433d177a5cSEmil Constantinescu if (adapt->ops->view) { 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ((*adapt->ops->view)(adapt,viewer)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 1473d177a5cSEmil Constantinescu } 1483d177a5cSEmil Constantinescu } 1493d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1503d177a5cSEmil Constantinescu } 1513d177a5cSEmil Constantinescu 1523d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt) 1533d177a5cSEmil Constantinescu { 1543d177a5cSEmil Constantinescu PetscFunctionBegin; 1553d177a5cSEmil Constantinescu if (!*adapt) PetscFunctionReturn(0); 1563d177a5cSEmil Constantinescu PetscValidHeaderSpecific(*adapt,TSGLLEADAPT_CLASSID,1); 157c793f718SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 1585f80ce2aSJacob Faibussowitsch if ((*adapt)->ops->destroy) CHKERRQ((*(*adapt)->ops->destroy)(*adapt)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHeaderDestroy(adapt)); 1603d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1613d177a5cSEmil Constantinescu } 1623d177a5cSEmil Constantinescu 1633d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSGLLEAdapt adapt) 1643d177a5cSEmil Constantinescu { 1653d177a5cSEmil Constantinescu char type[256] = TSGLLEADAPT_BOTH; 1663d177a5cSEmil Constantinescu PetscBool flg; 1673d177a5cSEmil Constantinescu 1683d177a5cSEmil Constantinescu PetscFunctionBegin; 1693d177a5cSEmil Constantinescu /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this 1703d177a5cSEmil Constantinescu * function can only be called from inside TSSetFromOptions_GLLE() */ 1715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"TSGLLE Adaptivity options")); 172*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSGLLEAdaptSetType",TSGLLEAdaptList, 173*b122ec5aSJacob Faibussowitsch ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg)); 174*b122ec5aSJacob Faibussowitsch if (flg || !((PetscObject)adapt)->type_name) CHKERRQ(TSGLLEAdaptSetType(adapt,type)); 1755f80ce2aSJacob Faibussowitsch if (adapt->ops->setfromoptions) CHKERRQ((*adapt->ops->setfromoptions)(PetscOptionsObject,adapt)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 1773d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1783d177a5cSEmil Constantinescu } 1793d177a5cSEmil Constantinescu 1803d177a5cSEmil 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) 1813d177a5cSEmil Constantinescu { 1823d177a5cSEmil Constantinescu PetscFunctionBegin; 1833d177a5cSEmil Constantinescu PetscValidHeaderSpecific(adapt,TSGLLEADAPT_CLASSID,1); 1843d177a5cSEmil Constantinescu PetscValidIntPointer(orders,3); 185dadcf809SJacob Faibussowitsch PetscValidRealPointer(errors,4); 186dadcf809SJacob Faibussowitsch PetscValidRealPointer(cost,5); 1873d177a5cSEmil Constantinescu PetscValidIntPointer(next_sc,9); 188dadcf809SJacob Faibussowitsch PetscValidRealPointer(next_h,10); 189064a246eSJacob Faibussowitsch PetscValidBoolPointer(finish,11); 1905f80ce2aSJacob Faibussowitsch CHKERRQ((*adapt->ops->choose)(adapt,n,orders,errors,cost,cur,h,tleft,next_sc,next_h,finish)); 1913d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1923d177a5cSEmil Constantinescu } 1933d177a5cSEmil Constantinescu 1943d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm,TSGLLEAdapt *inadapt) 1953d177a5cSEmil Constantinescu { 1963d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 1973d177a5cSEmil Constantinescu 1983d177a5cSEmil Constantinescu PetscFunctionBegin; 1993d177a5cSEmil Constantinescu *inadapt = NULL; 2005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHeaderCreate(adapt,TSGLLEADAPT_CLASSID,"TSGLLEAdapt","General Linear adaptivity","TS",comm,TSGLLEAdaptDestroy,TSGLLEAdaptView)); 2013d177a5cSEmil Constantinescu *inadapt = adapt; 2023d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2033d177a5cSEmil Constantinescu } 2043d177a5cSEmil Constantinescu 2053d177a5cSEmil Constantinescu /* 2060e3d61c9SBarry Smith Implementations 2073d177a5cSEmil Constantinescu */ 2083d177a5cSEmil Constantinescu 2093d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt) 2103d177a5cSEmil Constantinescu { 2113d177a5cSEmil Constantinescu PetscFunctionBegin; 2125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(adapt->data)); 2133d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2143d177a5cSEmil Constantinescu } 2153d177a5cSEmil Constantinescu 2163d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */ 2173d177a5cSEmil Constantinescu typedef struct { 2183d177a5cSEmil Constantinescu PetscInt scheme; 2193d177a5cSEmil Constantinescu PetscReal h; 2203d177a5cSEmil Constantinescu } TSGLLEAdapt_None; 2213d177a5cSEmil Constantinescu 2223d177a5cSEmil 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) 2233d177a5cSEmil Constantinescu { 2243d177a5cSEmil Constantinescu PetscFunctionBegin; 2253d177a5cSEmil Constantinescu *next_sc = cur; 2263d177a5cSEmil Constantinescu *next_h = h; 2273d177a5cSEmil Constantinescu if (*next_h > tleft) { 2283d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2293d177a5cSEmil Constantinescu *next_h = tleft; 2303d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2313d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2323d177a5cSEmil Constantinescu } 2333d177a5cSEmil Constantinescu 2343d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt) 2353d177a5cSEmil Constantinescu { 2363d177a5cSEmil Constantinescu TSGLLEAdapt_None *a; 2373d177a5cSEmil Constantinescu 2383d177a5cSEmil Constantinescu PetscFunctionBegin; 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(adapt,&a)); 2403d177a5cSEmil Constantinescu adapt->data = (void*)a; 2413d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_None; 2423d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2433d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2443d177a5cSEmil Constantinescu } 2453d177a5cSEmil Constantinescu 2463d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */ 2473d177a5cSEmil Constantinescu typedef struct { 2483d177a5cSEmil Constantinescu PetscReal desired_h; 2493d177a5cSEmil Constantinescu } TSGLLEAdapt_Size; 2503d177a5cSEmil Constantinescu 2513d177a5cSEmil 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) 2523d177a5cSEmil Constantinescu { 2533d177a5cSEmil Constantinescu TSGLLEAdapt_Size *sz = (TSGLLEAdapt_Size*)adapt->data; 2543d177a5cSEmil Constantinescu PetscReal dec = 0.2,inc = 5.0,safe = 0.9,optimal,last_desired_h; 2553d177a5cSEmil Constantinescu 2563d177a5cSEmil Constantinescu PetscFunctionBegin; 2573d177a5cSEmil Constantinescu *next_sc = cur; 2583d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[cur],(PetscReal)-1./(safe*orders[cur])); 2593d177a5cSEmil Constantinescu /* Step sizes oscillate when there is no smoothing. Here we use a geometric mean of the current step size and the 2603d177a5cSEmil Constantinescu * one that would have been taken (without smoothing) on the last step. */ 2613d177a5cSEmil Constantinescu last_desired_h = sz->desired_h; 2623d177a5cSEmil Constantinescu sz->desired_h = h*PetscMax(dec,PetscMin(inc,optimal)); /* Trim to [dec,inc] */ 2633d177a5cSEmil Constantinescu 2643d177a5cSEmil Constantinescu /* Normally only happens on the first step */ 2653d177a5cSEmil Constantinescu if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h); 2663d177a5cSEmil Constantinescu else *next_h = sz->desired_h; 2673d177a5cSEmil Constantinescu 2683d177a5cSEmil Constantinescu if (*next_h > tleft) { 2693d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2703d177a5cSEmil Constantinescu *next_h = tleft; 2713d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2723d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2733d177a5cSEmil Constantinescu } 2743d177a5cSEmil Constantinescu 2753d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt) 2763d177a5cSEmil Constantinescu { 2773d177a5cSEmil Constantinescu TSGLLEAdapt_Size *a; 2783d177a5cSEmil Constantinescu 2793d177a5cSEmil Constantinescu PetscFunctionBegin; 2805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(adapt,&a)); 2813d177a5cSEmil Constantinescu adapt->data = (void*)a; 2823d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Size; 2833d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2843d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2853d177a5cSEmil Constantinescu } 2863d177a5cSEmil Constantinescu 2873d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */ 2883d177a5cSEmil Constantinescu typedef struct { 2893d177a5cSEmil Constantinescu PetscInt count_at_order; 2903d177a5cSEmil Constantinescu PetscReal desired_h; 2913d177a5cSEmil Constantinescu } TSGLLEAdapt_Both; 2923d177a5cSEmil Constantinescu 2933d177a5cSEmil 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) 2943d177a5cSEmil Constantinescu { 2953d177a5cSEmil Constantinescu TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both*)adapt->data; 2963d177a5cSEmil Constantinescu PetscReal dec = 0.2,inc = 5.0,safe = 0.9; 2973d177a5cSEmil Constantinescu struct {PetscInt id; PetscReal h,eff;} best={-1,0,0},trial={-1,0,0},current={-1,0,0}; 2983d177a5cSEmil Constantinescu PetscInt i; 2993d177a5cSEmil Constantinescu 3003d177a5cSEmil Constantinescu PetscFunctionBegin; 3013d177a5cSEmil Constantinescu for (i=0; i<n; i++) { 3023d177a5cSEmil Constantinescu PetscReal optimal; 3033d177a5cSEmil Constantinescu trial.id = i; 3043d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[i],(PetscReal)-1./(safe*orders[i])); 3053d177a5cSEmil Constantinescu trial.h = h*optimal; 3063d177a5cSEmil Constantinescu trial.eff = trial.h/cost[i]; 3075f80ce2aSJacob Faibussowitsch if (trial.eff > best.eff) CHKERRQ(PetscArraycpy(&best,&trial,1)); 3085f80ce2aSJacob Faibussowitsch if (i == cur) CHKERRQ(PetscArraycpy(¤t,&trial,1)); 3093d177a5cSEmil Constantinescu } 3103d177a5cSEmil Constantinescu /* Only switch orders if the scheme offers significant benefits over the current one. 3113d177a5cSEmil Constantinescu When the scheme is not changing, only change step size if it offers significant benefits. */ 3123d177a5cSEmil Constantinescu if (best.eff < 1.2*current.eff || both->count_at_order < orders[cur]+2) { 3133d177a5cSEmil Constantinescu PetscReal last_desired_h; 3143d177a5cSEmil Constantinescu *next_sc = current.id; 3153d177a5cSEmil Constantinescu last_desired_h = both->desired_h; 3163d177a5cSEmil Constantinescu both->desired_h = PetscMax(h*dec,PetscMin(h*inc,current.h)); 3173d177a5cSEmil Constantinescu *next_h = (both->count_at_order > 0) 3183d177a5cSEmil Constantinescu ? PetscSqrtReal(last_desired_h * both->desired_h) 3193d177a5cSEmil Constantinescu : both->desired_h; 3203d177a5cSEmil Constantinescu both->count_at_order++; 3213d177a5cSEmil Constantinescu } else { 3223d177a5cSEmil Constantinescu PetscReal rat = cost[best.id]/cost[cur]; 3233d177a5cSEmil Constantinescu *next_sc = best.id; 3243d177a5cSEmil Constantinescu *next_h = PetscMax(h*rat*dec,PetscMin(h*rat*inc,best.h)); 3253d177a5cSEmil Constantinescu both->count_at_order = 0; 3263d177a5cSEmil Constantinescu both->desired_h = best.h; 3273d177a5cSEmil Constantinescu } 3283d177a5cSEmil Constantinescu 3293d177a5cSEmil Constantinescu if (*next_h > tleft) { 3303d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 3313d177a5cSEmil Constantinescu *next_h = tleft; 3323d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 3333d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3343d177a5cSEmil Constantinescu } 3353d177a5cSEmil Constantinescu 3363d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt) 3373d177a5cSEmil Constantinescu { 3383d177a5cSEmil Constantinescu TSGLLEAdapt_Both *a; 3393d177a5cSEmil Constantinescu 3403d177a5cSEmil Constantinescu PetscFunctionBegin; 3415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(adapt,&a)); 3423d177a5cSEmil Constantinescu adapt->data = (void*)a; 3433d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Both; 3443d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 3453d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3463d177a5cSEmil Constantinescu } 347