1*3d177a5cSEmil Constantinescu 2*3d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I "petscts.h" I*/ 3*3d177a5cSEmil Constantinescu 4*3d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAdaptList; 5*3d177a5cSEmil Constantinescu static PetscBool TSGLLEAdaptPackageInitialized; 6*3d177a5cSEmil Constantinescu static PetscBool TSGLLEAdaptRegisterAllCalled; 7*3d177a5cSEmil Constantinescu static PetscClassId TSGLLEADAPT_CLASSID; 8*3d177a5cSEmil Constantinescu 9*3d177a5cSEmil Constantinescu struct _TSGLLEAdaptOps { 10*3d177a5cSEmil Constantinescu PetscErrorCode (*choose)(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool*); 11*3d177a5cSEmil Constantinescu PetscErrorCode (*destroy)(TSGLLEAdapt); 12*3d177a5cSEmil Constantinescu PetscErrorCode (*view)(TSGLLEAdapt,PetscViewer); 13*3d177a5cSEmil Constantinescu PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSGLLEAdapt); 14*3d177a5cSEmil Constantinescu }; 15*3d177a5cSEmil Constantinescu 16*3d177a5cSEmil Constantinescu struct _p_TSGLLEAdapt { 17*3d177a5cSEmil Constantinescu PETSCHEADER(struct _TSGLLEAdaptOps); 18*3d177a5cSEmil Constantinescu void *data; 19*3d177a5cSEmil Constantinescu }; 20*3d177a5cSEmil Constantinescu 21*3d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt); 22*3d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt); 23*3d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt); 24*3d177a5cSEmil Constantinescu 25*3d177a5cSEmil Constantinescu #undef __FUNCT__ 26*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptRegister" 27*3d177a5cSEmil Constantinescu /*@C 28*3d177a5cSEmil Constantinescu TSGLLEAdaptRegister - adds a TSGLLEAdapt implementation 29*3d177a5cSEmil Constantinescu 30*3d177a5cSEmil Constantinescu Not Collective 31*3d177a5cSEmil Constantinescu 32*3d177a5cSEmil Constantinescu Input Parameters: 33*3d177a5cSEmil Constantinescu + name_scheme - name of user-defined adaptivity scheme 34*3d177a5cSEmil Constantinescu - routine_create - routine to create method context 35*3d177a5cSEmil Constantinescu 36*3d177a5cSEmil Constantinescu Notes: 37*3d177a5cSEmil Constantinescu TSGLLEAdaptRegister() may be called multiple times to add several user-defined families. 38*3d177a5cSEmil Constantinescu 39*3d177a5cSEmil Constantinescu Sample usage: 40*3d177a5cSEmil Constantinescu .vb 41*3d177a5cSEmil Constantinescu TSGLLEAdaptRegister("my_scheme",MySchemeCreate); 42*3d177a5cSEmil Constantinescu .ve 43*3d177a5cSEmil Constantinescu 44*3d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 45*3d177a5cSEmil Constantinescu $ TSGLLEAdaptSetType(ts,"my_scheme") 46*3d177a5cSEmil Constantinescu or at runtime via the option 47*3d177a5cSEmil Constantinescu $ -ts_adapt_type my_scheme 48*3d177a5cSEmil Constantinescu 49*3d177a5cSEmil Constantinescu Level: advanced 50*3d177a5cSEmil Constantinescu 51*3d177a5cSEmil Constantinescu .keywords: TSGLLEAdapt, register 52*3d177a5cSEmil Constantinescu 53*3d177a5cSEmil Constantinescu .seealso: TSGLLEAdaptRegisterAll() 54*3d177a5cSEmil Constantinescu @*/ 55*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptRegister(const char sname[],PetscErrorCode (*function)(TSGLLEAdapt)) 56*3d177a5cSEmil Constantinescu { 57*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 58*3d177a5cSEmil Constantinescu 59*3d177a5cSEmil Constantinescu PetscFunctionBegin; 60*3d177a5cSEmil Constantinescu ierr = PetscFunctionListAdd(&TSGLLEAdaptList,sname,function);CHKERRQ(ierr); 61*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 62*3d177a5cSEmil Constantinescu } 63*3d177a5cSEmil Constantinescu 64*3d177a5cSEmil Constantinescu #undef __FUNCT__ 65*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptRegisterAll" 66*3d177a5cSEmil Constantinescu /*@C 67*3d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in TSGLLEAdapt 68*3d177a5cSEmil Constantinescu 69*3d177a5cSEmil Constantinescu Not Collective 70*3d177a5cSEmil Constantinescu 71*3d177a5cSEmil Constantinescu Level: advanced 72*3d177a5cSEmil Constantinescu 73*3d177a5cSEmil Constantinescu .keywords: TSGLLEAdapt, register, all 74*3d177a5cSEmil Constantinescu 75*3d177a5cSEmil Constantinescu .seealso: TSGLLEAdaptRegisterDestroy() 76*3d177a5cSEmil Constantinescu @*/ 77*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptRegisterAll(void) 78*3d177a5cSEmil Constantinescu { 79*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 80*3d177a5cSEmil Constantinescu 81*3d177a5cSEmil Constantinescu PetscFunctionBegin; 82*3d177a5cSEmil Constantinescu if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(0); 83*3d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_TRUE; 84*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegister(TSGLLEADAPT_NONE,TSGLLEAdaptCreate_None);CHKERRQ(ierr); 85*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegister(TSGLLEADAPT_SIZE,TSGLLEAdaptCreate_Size);CHKERRQ(ierr); 86*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegister(TSGLLEADAPT_BOTH,TSGLLEAdaptCreate_Both);CHKERRQ(ierr); 87*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 88*3d177a5cSEmil Constantinescu } 89*3d177a5cSEmil Constantinescu 90*3d177a5cSEmil Constantinescu #undef __FUNCT__ 91*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptFinalizePackage" 92*3d177a5cSEmil Constantinescu /*@C 93*3d177a5cSEmil Constantinescu TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 94*3d177a5cSEmil Constantinescu called from PetscFinalize(). 95*3d177a5cSEmil Constantinescu 96*3d177a5cSEmil Constantinescu Level: developer 97*3d177a5cSEmil Constantinescu 98*3d177a5cSEmil Constantinescu .keywords: Petsc, destroy, package 99*3d177a5cSEmil Constantinescu .seealso: PetscFinalize() 100*3d177a5cSEmil Constantinescu @*/ 101*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptFinalizePackage(void) 102*3d177a5cSEmil Constantinescu { 103*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 104*3d177a5cSEmil Constantinescu 105*3d177a5cSEmil Constantinescu PetscFunctionBegin; 106*3d177a5cSEmil Constantinescu ierr = PetscFunctionListDestroy(&TSGLLEAdaptList);CHKERRQ(ierr); 107*3d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_FALSE; 108*3d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_FALSE; 109*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 110*3d177a5cSEmil Constantinescu } 111*3d177a5cSEmil Constantinescu 112*3d177a5cSEmil Constantinescu #undef __FUNCT__ 113*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptInitializePackage" 114*3d177a5cSEmil Constantinescu /*@C 115*3d177a5cSEmil Constantinescu TSGLLEAdaptInitializePackage - This function initializes everything in the TSGLLEAdapt package. It is 116*3d177a5cSEmil Constantinescu called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 117*3d177a5cSEmil Constantinescu TSCreate_GLLE() when using static libraries. 118*3d177a5cSEmil Constantinescu 119*3d177a5cSEmil Constantinescu Level: developer 120*3d177a5cSEmil Constantinescu 121*3d177a5cSEmil Constantinescu .keywords: TSGLLEAdapt, initialize, package 122*3d177a5cSEmil Constantinescu .seealso: PetscInitialize() 123*3d177a5cSEmil Constantinescu @*/ 124*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptInitializePackage(void) 125*3d177a5cSEmil Constantinescu { 126*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 127*3d177a5cSEmil Constantinescu 128*3d177a5cSEmil Constantinescu PetscFunctionBegin; 129*3d177a5cSEmil Constantinescu if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(0); 130*3d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_TRUE; 131*3d177a5cSEmil Constantinescu ierr = PetscClassIdRegister("TSGLLEAdapt",&TSGLLEADAPT_CLASSID);CHKERRQ(ierr); 132*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptRegisterAll();CHKERRQ(ierr); 133*3d177a5cSEmil Constantinescu ierr = PetscRegisterFinalize(TSGLLEAdaptFinalizePackage);CHKERRQ(ierr); 134*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 135*3d177a5cSEmil Constantinescu } 136*3d177a5cSEmil Constantinescu 137*3d177a5cSEmil Constantinescu #undef __FUNCT__ 138*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptSetType" 139*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt,TSGLLEAdaptType type) 140*3d177a5cSEmil Constantinescu { 141*3d177a5cSEmil Constantinescu PetscErrorCode ierr,(*r)(TSGLLEAdapt); 142*3d177a5cSEmil Constantinescu 143*3d177a5cSEmil Constantinescu PetscFunctionBegin; 144*3d177a5cSEmil Constantinescu ierr = PetscFunctionListFind(TSGLLEAdaptList,type,&r);CHKERRQ(ierr); 145*3d177a5cSEmil Constantinescu if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAdapt type \"%s\" given",type); 146*3d177a5cSEmil Constantinescu if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 147*3d177a5cSEmil Constantinescu ierr = (*r)(adapt);CHKERRQ(ierr); 148*3d177a5cSEmil Constantinescu ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 149*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 150*3d177a5cSEmil Constantinescu } 151*3d177a5cSEmil Constantinescu 152*3d177a5cSEmil Constantinescu #undef __FUNCT__ 153*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptSetOptionsPrefix" 154*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt,const char prefix[]) 155*3d177a5cSEmil Constantinescu { 156*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 157*3d177a5cSEmil Constantinescu 158*3d177a5cSEmil Constantinescu PetscFunctionBegin; 159*3d177a5cSEmil Constantinescu ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 160*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 161*3d177a5cSEmil Constantinescu } 162*3d177a5cSEmil Constantinescu 163*3d177a5cSEmil Constantinescu #undef __FUNCT__ 164*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptView" 165*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt,PetscViewer viewer) 166*3d177a5cSEmil Constantinescu { 167*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 168*3d177a5cSEmil Constantinescu PetscBool iascii; 169*3d177a5cSEmil Constantinescu 170*3d177a5cSEmil Constantinescu PetscFunctionBegin; 171*3d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 172*3d177a5cSEmil Constantinescu if (iascii) { 173*3d177a5cSEmil Constantinescu ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 174*3d177a5cSEmil Constantinescu if (adapt->ops->view) { 175*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 176*3d177a5cSEmil Constantinescu ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 177*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 178*3d177a5cSEmil Constantinescu } 179*3d177a5cSEmil Constantinescu } 180*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 181*3d177a5cSEmil Constantinescu } 182*3d177a5cSEmil Constantinescu 183*3d177a5cSEmil Constantinescu #undef __FUNCT__ 184*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptDestroy" 185*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt) 186*3d177a5cSEmil Constantinescu { 187*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 188*3d177a5cSEmil Constantinescu 189*3d177a5cSEmil Constantinescu PetscFunctionBegin; 190*3d177a5cSEmil Constantinescu if (!*adapt) PetscFunctionReturn(0); 191*3d177a5cSEmil Constantinescu PetscValidHeaderSpecific(*adapt,TSGLLEADAPT_CLASSID,1); 192*3d177a5cSEmil Constantinescu if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);} 193*3d177a5cSEmil Constantinescu if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 194*3d177a5cSEmil Constantinescu ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 195*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 196*3d177a5cSEmil Constantinescu } 197*3d177a5cSEmil Constantinescu 198*3d177a5cSEmil Constantinescu #undef __FUNCT__ 199*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptSetFromOptions" 200*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSGLLEAdapt adapt) 201*3d177a5cSEmil Constantinescu { 202*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 203*3d177a5cSEmil Constantinescu char type[256] = TSGLLEADAPT_BOTH; 204*3d177a5cSEmil Constantinescu PetscBool flg; 205*3d177a5cSEmil Constantinescu 206*3d177a5cSEmil Constantinescu PetscFunctionBegin; 207*3d177a5cSEmil Constantinescu /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this 208*3d177a5cSEmil Constantinescu * function can only be called from inside TSSetFromOptions_GLLE() */ 209*3d177a5cSEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"TSGLLE Adaptivity options");CHKERRQ(ierr); 210*3d177a5cSEmil Constantinescu ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSGLLEAdaptSetType",TSGLLEAdaptList, 211*3d177a5cSEmil Constantinescu ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 212*3d177a5cSEmil Constantinescu if (flg || !((PetscObject)adapt)->type_name) { 213*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptSetType(adapt,type);CHKERRQ(ierr); 214*3d177a5cSEmil Constantinescu } 215*3d177a5cSEmil Constantinescu if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 216*3d177a5cSEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 217*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 218*3d177a5cSEmil Constantinescu } 219*3d177a5cSEmil Constantinescu 220*3d177a5cSEmil Constantinescu #undef __FUNCT__ 221*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptChoose" 222*3d177a5cSEmil 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) 223*3d177a5cSEmil Constantinescu { 224*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 225*3d177a5cSEmil Constantinescu 226*3d177a5cSEmil Constantinescu PetscFunctionBegin; 227*3d177a5cSEmil Constantinescu PetscValidHeaderSpecific(adapt,TSGLLEADAPT_CLASSID,1); 228*3d177a5cSEmil Constantinescu PetscValidIntPointer(orders,3); 229*3d177a5cSEmil Constantinescu PetscValidPointer(errors,4); 230*3d177a5cSEmil Constantinescu PetscValidPointer(cost,5); 231*3d177a5cSEmil Constantinescu PetscValidIntPointer(next_sc,9); 232*3d177a5cSEmil Constantinescu PetscValidPointer(next_h,10); 233*3d177a5cSEmil Constantinescu PetscValidIntPointer(finish,11); 234*3d177a5cSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,n,orders,errors,cost,cur,h,tleft,next_sc,next_h,finish);CHKERRQ(ierr); 235*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 236*3d177a5cSEmil Constantinescu } 237*3d177a5cSEmil Constantinescu 238*3d177a5cSEmil Constantinescu #undef __FUNCT__ 239*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptCreate" 240*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm,TSGLLEAdapt *inadapt) 241*3d177a5cSEmil Constantinescu { 242*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 243*3d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 244*3d177a5cSEmil Constantinescu 245*3d177a5cSEmil Constantinescu PetscFunctionBegin; 246*3d177a5cSEmil Constantinescu *inadapt = NULL; 247*3d177a5cSEmil Constantinescu ierr = PetscHeaderCreate(adapt,TSGLLEADAPT_CLASSID,"TSGLLEAdapt","General Linear adaptivity","TS",comm,TSGLLEAdaptDestroy,TSGLLEAdaptView);CHKERRQ(ierr); 248*3d177a5cSEmil Constantinescu *inadapt = adapt; 249*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 250*3d177a5cSEmil Constantinescu } 251*3d177a5cSEmil Constantinescu 252*3d177a5cSEmil Constantinescu 253*3d177a5cSEmil Constantinescu /* 254*3d177a5cSEmil Constantinescu * Implementations 255*3d177a5cSEmil Constantinescu */ 256*3d177a5cSEmil Constantinescu 257*3d177a5cSEmil Constantinescu #undef __FUNCT__ 258*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptDestroy_JustFree" 259*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt) 260*3d177a5cSEmil Constantinescu { 261*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 262*3d177a5cSEmil Constantinescu 263*3d177a5cSEmil Constantinescu PetscFunctionBegin; 264*3d177a5cSEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 265*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 266*3d177a5cSEmil Constantinescu } 267*3d177a5cSEmil Constantinescu 268*3d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */ 269*3d177a5cSEmil Constantinescu typedef struct { 270*3d177a5cSEmil Constantinescu PetscInt scheme; 271*3d177a5cSEmil Constantinescu PetscReal h; 272*3d177a5cSEmil Constantinescu } TSGLLEAdapt_None; 273*3d177a5cSEmil Constantinescu 274*3d177a5cSEmil Constantinescu #undef __FUNCT__ 275*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptChoose_None" 276*3d177a5cSEmil 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) 277*3d177a5cSEmil Constantinescu { 278*3d177a5cSEmil Constantinescu 279*3d177a5cSEmil Constantinescu PetscFunctionBegin; 280*3d177a5cSEmil Constantinescu *next_sc = cur; 281*3d177a5cSEmil Constantinescu *next_h = h; 282*3d177a5cSEmil Constantinescu if (*next_h > tleft) { 283*3d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 284*3d177a5cSEmil Constantinescu *next_h = tleft; 285*3d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 286*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 287*3d177a5cSEmil Constantinescu } 288*3d177a5cSEmil Constantinescu 289*3d177a5cSEmil Constantinescu #undef __FUNCT__ 290*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptCreate_None" 291*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt) 292*3d177a5cSEmil Constantinescu { 293*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 294*3d177a5cSEmil Constantinescu TSGLLEAdapt_None *a; 295*3d177a5cSEmil Constantinescu 296*3d177a5cSEmil Constantinescu PetscFunctionBegin; 297*3d177a5cSEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 298*3d177a5cSEmil Constantinescu adapt->data = (void*)a; 299*3d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_None; 300*3d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 301*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 302*3d177a5cSEmil Constantinescu } 303*3d177a5cSEmil Constantinescu 304*3d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */ 305*3d177a5cSEmil Constantinescu typedef struct { 306*3d177a5cSEmil Constantinescu PetscReal desired_h; 307*3d177a5cSEmil Constantinescu } TSGLLEAdapt_Size; 308*3d177a5cSEmil Constantinescu 309*3d177a5cSEmil Constantinescu 310*3d177a5cSEmil Constantinescu #undef __FUNCT__ 311*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptChoose_Size" 312*3d177a5cSEmil 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) 313*3d177a5cSEmil Constantinescu { 314*3d177a5cSEmil Constantinescu TSGLLEAdapt_Size *sz = (TSGLLEAdapt_Size*)adapt->data; 315*3d177a5cSEmil Constantinescu PetscReal dec = 0.2,inc = 5.0,safe = 0.9,optimal,last_desired_h; 316*3d177a5cSEmil Constantinescu 317*3d177a5cSEmil Constantinescu PetscFunctionBegin; 318*3d177a5cSEmil Constantinescu *next_sc = cur; 319*3d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[cur],(PetscReal)-1./(safe*orders[cur])); 320*3d177a5cSEmil Constantinescu /* Step sizes oscillate when there is no smoothing. Here we use a geometric mean of the current step size and the 321*3d177a5cSEmil Constantinescu * one that would have been taken (without smoothing) on the last step. */ 322*3d177a5cSEmil Constantinescu last_desired_h = sz->desired_h; 323*3d177a5cSEmil Constantinescu sz->desired_h = h*PetscMax(dec,PetscMin(inc,optimal)); /* Trim to [dec,inc] */ 324*3d177a5cSEmil Constantinescu 325*3d177a5cSEmil Constantinescu /* Normally only happens on the first step */ 326*3d177a5cSEmil Constantinescu if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h); 327*3d177a5cSEmil Constantinescu else *next_h = sz->desired_h; 328*3d177a5cSEmil Constantinescu 329*3d177a5cSEmil Constantinescu if (*next_h > tleft) { 330*3d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 331*3d177a5cSEmil Constantinescu *next_h = tleft; 332*3d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 333*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 334*3d177a5cSEmil Constantinescu } 335*3d177a5cSEmil Constantinescu 336*3d177a5cSEmil Constantinescu #undef __FUNCT__ 337*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptCreate_Size" 338*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt) 339*3d177a5cSEmil Constantinescu { 340*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 341*3d177a5cSEmil Constantinescu TSGLLEAdapt_Size *a; 342*3d177a5cSEmil Constantinescu 343*3d177a5cSEmil Constantinescu PetscFunctionBegin; 344*3d177a5cSEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 345*3d177a5cSEmil Constantinescu adapt->data = (void*)a; 346*3d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Size; 347*3d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 348*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 349*3d177a5cSEmil Constantinescu } 350*3d177a5cSEmil Constantinescu 351*3d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */ 352*3d177a5cSEmil Constantinescu typedef struct { 353*3d177a5cSEmil Constantinescu PetscInt count_at_order; 354*3d177a5cSEmil Constantinescu PetscReal desired_h; 355*3d177a5cSEmil Constantinescu } TSGLLEAdapt_Both; 356*3d177a5cSEmil Constantinescu 357*3d177a5cSEmil Constantinescu 358*3d177a5cSEmil Constantinescu #undef __FUNCT__ 359*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptChoose_Both" 360*3d177a5cSEmil 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) 361*3d177a5cSEmil Constantinescu { 362*3d177a5cSEmil Constantinescu TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both*)adapt->data; 363*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 364*3d177a5cSEmil Constantinescu PetscReal dec = 0.2,inc = 5.0,safe = 0.9; 365*3d177a5cSEmil Constantinescu struct {PetscInt id; PetscReal h,eff;} best={-1,0,0},trial={-1,0,0},current={-1,0,0}; 366*3d177a5cSEmil Constantinescu PetscInt i; 367*3d177a5cSEmil Constantinescu 368*3d177a5cSEmil Constantinescu PetscFunctionBegin; 369*3d177a5cSEmil Constantinescu for (i=0; i<n; i++) { 370*3d177a5cSEmil Constantinescu PetscReal optimal; 371*3d177a5cSEmil Constantinescu trial.id = i; 372*3d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[i],(PetscReal)-1./(safe*orders[i])); 373*3d177a5cSEmil Constantinescu trial.h = h*optimal; 374*3d177a5cSEmil Constantinescu trial.eff = trial.h/cost[i]; 375*3d177a5cSEmil Constantinescu if (trial.eff > best.eff) {ierr = PetscMemcpy(&best,&trial,sizeof(trial));CHKERRQ(ierr);} 376*3d177a5cSEmil Constantinescu if (i == cur) {ierr = PetscMemcpy(¤t,&trial,sizeof(trial));CHKERRQ(ierr);} 377*3d177a5cSEmil Constantinescu } 378*3d177a5cSEmil Constantinescu /* Only switch orders if the scheme offers significant benefits over the current one. 379*3d177a5cSEmil Constantinescu When the scheme is not changing, only change step size if it offers significant benefits. */ 380*3d177a5cSEmil Constantinescu if (best.eff < 1.2*current.eff || both->count_at_order < orders[cur]+2) { 381*3d177a5cSEmil Constantinescu PetscReal last_desired_h; 382*3d177a5cSEmil Constantinescu *next_sc = current.id; 383*3d177a5cSEmil Constantinescu last_desired_h = both->desired_h; 384*3d177a5cSEmil Constantinescu both->desired_h = PetscMax(h*dec,PetscMin(h*inc,current.h)); 385*3d177a5cSEmil Constantinescu *next_h = (both->count_at_order > 0) 386*3d177a5cSEmil Constantinescu ? PetscSqrtReal(last_desired_h * both->desired_h) 387*3d177a5cSEmil Constantinescu : both->desired_h; 388*3d177a5cSEmil Constantinescu both->count_at_order++; 389*3d177a5cSEmil Constantinescu } else { 390*3d177a5cSEmil Constantinescu PetscReal rat = cost[best.id]/cost[cur]; 391*3d177a5cSEmil Constantinescu *next_sc = best.id; 392*3d177a5cSEmil Constantinescu *next_h = PetscMax(h*rat*dec,PetscMin(h*rat*inc,best.h)); 393*3d177a5cSEmil Constantinescu both->count_at_order = 0; 394*3d177a5cSEmil Constantinescu both->desired_h = best.h; 395*3d177a5cSEmil Constantinescu } 396*3d177a5cSEmil Constantinescu 397*3d177a5cSEmil Constantinescu if (*next_h > tleft) { 398*3d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 399*3d177a5cSEmil Constantinescu *next_h = tleft; 400*3d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 401*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 402*3d177a5cSEmil Constantinescu } 403*3d177a5cSEmil Constantinescu 404*3d177a5cSEmil Constantinescu #undef __FUNCT__ 405*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAdaptCreate_Both" 406*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt) 407*3d177a5cSEmil Constantinescu { 408*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 409*3d177a5cSEmil Constantinescu TSGLLEAdapt_Both *a; 410*3d177a5cSEmil Constantinescu 411*3d177a5cSEmil Constantinescu PetscFunctionBegin; 412*3d177a5cSEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 413*3d177a5cSEmil Constantinescu adapt->data = (void*)a; 414*3d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Both; 415*3d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 416*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 417*3d177a5cSEmil Constantinescu } 418