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