Coverage for pySDC/implementations/convergence_controller_classes/adaptivity.py: 95%

257 statements  

« prev     ^ index     » next       coverage.py v7.5.0, created at 2024-04-29 09:02 +0000

1import numpy as np 

2from pySDC.core.ConvergenceController import ConvergenceController, Status 

3from pySDC.implementations.convergence_controller_classes.step_size_limiter import ( 

4 StepSizeLimiter, 

5) 

6 

7 

8class AdaptivityBase(ConvergenceController): 

9 """ 

10 Abstract base class for convergence controllers that implement adaptivity based on arbitrary local error estimates 

11 and update rules. 

12 """ 

13 

14 def setup(self, controller, params, description, **kwargs): 

15 """ 

16 Define default parameters here. 

17 

18 Default parameters are: 

19 - control_order (int): The order relative to other convergence controllers 

20 - beta (float): The safety factor 

21 

22 Args: 

23 controller (pySDC.Controller): The controller 

24 params (dict): The params passed for this specific convergence controller 

25 description (dict): The description object used to instantiate the controller 

26 

27 Returns: 

28 (dict): The updated params dictionary 

29 """ 

30 defaults = { 

31 "control_order": -50, 

32 "beta": 0.9, 

33 } 

34 from pySDC.implementations.hooks.log_step_size import LogStepSize 

35 

36 controller.add_hook(LogStepSize) 

37 

38 from pySDC.implementations.convergence_controller_classes.check_convergence import CheckConvergence 

39 

40 self.communicate_convergence = CheckConvergence.communicate_convergence 

41 

42 return {**defaults, **super().setup(controller, params, description, **kwargs)} 

43 

44 def dependencies(self, controller, description, **kwargs): 

45 """ 

46 Load step size limiters here, if they are desired. 

47 

48 Args: 

49 controller (pySDC.Controller): The controller 

50 description (dict): The description object used to instantiate the controller 

51 

52 Returns: 

53 None 

54 """ 

55 step_limiter_keys = ['dt_min', 'dt_max', 'dt_slope_min', 'dt_slope_max'] 

56 available_keys = [me for me in step_limiter_keys if me in self.params.__dict__.keys()] 

57 

58 if len(available_keys) > 0: 

59 step_limiter_params = {key: self.params.__dict__[key] for key in available_keys} 

60 controller.add_convergence_controller(StepSizeLimiter, params=step_limiter_params, description=description) 

61 

62 if self.params.useMPI: 

63 self.prepare_MPI_logical_operations() 

64 

65 return None 

66 

67 def get_new_step_size(self, controller, S, **kwargs): 

68 """ 

69 Determine a step size for the next step from an estimate of the local error of the current step. 

70 

71 Args: 

72 controller (pySDC.Controller): The controller 

73 S (pySDC.Step): The current step 

74 

75 Returns: 

76 None 

77 """ 

78 raise NotImplementedError("Please implement a rule for updating the step size!") 

79 

80 def compute_optimal_step_size(self, beta, dt, e_tol, e_est, order): 

81 """ 

82 Compute the optimal step size for the current step based on the order of the scheme. 

83 This function can be called from `get_new_step_size` for various implementations of adaptivity, but notably not 

84 all! We require to know the order of the error estimate and if we do adaptivity based on the residual, for 

85 instance, we do not know that and we can't use this function. 

86 

87 Args: 

88 beta (float): Safety factor 

89 dt (float): Current step size 

90 e_tol (float): The desired tolerance 

91 e_est (float): The estimated local error 

92 order (int): The order of the local error estimate 

93 

94 Returns: 

95 float: The optimal step size 

96 """ 

97 return beta * dt * (e_tol / e_est) ** (1.0 / order) 

98 

99 def get_local_error_estimate(self, controller, S, **kwargs): 

100 """ 

101 Get the local error estimate for updating the step size. 

102 It does not have to be an error estimate, but could be the residual or something else. 

103 

104 Args: 

105 controller (pySDC.Controller): The controller 

106 S (pySDC.Step): The current step 

107 

108 Returns: 

109 float: The error estimate 

110 """ 

111 raise NotImplementedError("Please implement a way to get the local error") 

112 

113 def determine_restart(self, controller, S, **kwargs): 

114 """ 

115 Check if the step wants to be restarted by comparing the estimate of the local error to a preset tolerance 

116 

117 Args: 

118 controller (pySDC.Controller): The controller 

119 S (pySDC.Step): The current step 

120 

121 Returns: 

122 None 

123 """ 

124 if S.status.iter >= S.params.maxiter: 

125 e_est = self.get_local_error_estimate(controller, S) 

126 if e_est >= self.params.e_tol: 

127 # see if we try to avoid restarts 

128 if self.params.get('avoid_restarts'): 

129 more_iter_needed = max([L.status.iter_to_convergence for L in S.levels]) 

130 k_final = S.status.iter + more_iter_needed 

131 rho = max([L.status.contraction_factor for L in S.levels]) 

132 coll_order = S.levels[0].sweep.coll.order 

133 

134 if rho > 1: 

135 S.status.restart = True 

136 self.log(f"Convergence factor = {rho:.2e} > 1 -> restarting", S) 

137 elif k_final > 2 * S.params.maxiter: 

138 S.status.restart = True 

139 self.log( 

140 f"{more_iter_needed} more iterations needed for convergence -> restart is more efficient", S 

141 ) 

142 elif k_final > coll_order: 

143 S.status.restart = True 

144 self.log( 

145 f"{more_iter_needed} more iterations needed for convergence -> restart because collocation problem would be over resolved", 

146 S, 

147 ) 

148 else: 

149 S.status.force_continue = True 

150 self.log(f"{more_iter_needed} more iterations needed for convergence -> no restart", S) 

151 else: 

152 S.status.restart = True 

153 self.log(f"Restarting: e={e_est:.2e} >= e_tol={self.params.e_tol:.2e}", S) 

154 

155 return None 

156 

157 

158class AdaptivityForConvergedCollocationProblems(AdaptivityBase): 

159 def dependencies(self, controller, description, **kwargs): 

160 """ 

161 Load interpolation between restarts. 

162 

163 Args: 

164 controller (pySDC.Controller): The controller 

165 description (dict): The description object used to instantiate the controller 

166 

167 Returns: 

168 None 

169 """ 

170 super().dependencies(controller, description, **kwargs) 

171 

172 if self.params.interpolate_between_restarts: 

173 from pySDC.implementations.convergence_controller_classes.interpolate_between_restarts import ( 

174 InterpolateBetweenRestarts, 

175 ) 

176 

177 controller.add_convergence_controller( 

178 InterpolateBetweenRestarts, 

179 description=description, 

180 params={}, 

181 ) 

182 if self.params.interpolate_between_restarts: 

183 self.interpolator = controller.convergence_controllers[-1] 

184 return None 

185 

186 def get_convergence(self, controller, S, **kwargs): 

187 raise NotImplementedError("Please implement a way to check if the collocation problem is converged!") 

188 

189 def setup(self, controller, params, description, **kwargs): 

190 """ 

191 Add a default value for control order to the parameters. 

192 

193 Args: 

194 controller (pySDC.Controller): The controller 

195 params (dict): Parameters for the convergence controller 

196 description (dict): The description object used to instantiate the controller 

197 

198 Returns: 

199 dict: Updated parameters 

200 """ 

201 defaults = { 

202 'restol_rel': None, 

203 'e_tol_rel': None, 

204 'restart_at_maxiter': True, 

205 'restol_min': 1e-12, 

206 'restol_max': 1e-5, 

207 'factor_if_not_converged': 4.0, 

208 'residual_max_tol': 1e9, 

209 'maxiter': description['sweeper_params'].get('maxiter', 99), 

210 'interpolate_between_restarts': True, 

211 'abort_at_growing_residual': True, 

212 **super().setup(controller, params, description, **kwargs), 

213 } 

214 if defaults['restol_rel']: 

215 description['level_params']['restol'] = min( 

216 [max([defaults['restol_rel'] * defaults['e_tol'], defaults['restol_min']]), defaults['restol_max']] 

217 ) 

218 elif defaults['e_tol_rel']: 

219 description['level_params']['e_tol'] = min([max([defaults['e_tol_rel'] * defaults['e_tol'], 1e-10]), 1e-5]) 

220 

221 if defaults['restart_at_maxiter']: 

222 defaults['maxiter'] = description['step_params'].get('maxiter', 99) 

223 

224 self.res_last_iter = np.inf 

225 

226 return defaults 

227 

228 def determine_restart(self, controller, S, **kwargs): 

229 if self.get_convergence(controller, S, **kwargs): 

230 self.res_last_iter = np.inf 

231 

232 if self.params.restart_at_maxiter and S.levels[0].status.residual > S.levels[0].params.restol: 

233 self.trigger_restart_upon_nonconvergence(S) 

234 elif self.get_local_error_estimate(controller, S, **kwargs) > self.params.e_tol: 

235 S.status.restart = True 

236 elif ( 

237 S.status.time_size == 1 

238 and self.res_last_iter < S.levels[0].status.residual 

239 and S.status.iter > 0 

240 and self.params.abort_at_growing_residual 

241 ): 

242 self.trigger_restart_upon_nonconvergence(S) 

243 elif S.levels[0].status.residual > self.params.residual_max_tol: 

244 self.trigger_restart_upon_nonconvergence(S) 

245 

246 if self.params.useMPI: 

247 self.communicate_convergence(self, controller, S, **kwargs) 

248 

249 self.res_last_iter = S.levels[0].status.residual * 1.0 

250 

251 def trigger_restart_upon_nonconvergence(self, S): 

252 S.status.restart = True 

253 S.status.force_done = True 

254 for L in S.levels: 

255 L.status.dt_new = L.params.dt / self.params.factor_if_not_converged 

256 self.log( 

257 f'Collocation problem not converged. Reducing step size to {L.status.dt_new:.2e}', 

258 S, 

259 ) 

260 if self.params.interpolate_between_restarts: 

261 self.interpolator.status.skip_interpolation = True 

262 

263 

264class Adaptivity(AdaptivityBase): 

265 """ 

266 Class to compute time step size adaptively based on embedded error estimate. 

267 

268 We have a version working in non-MPI pipelined SDC, but Adaptivity requires you to know the order of the scheme, 

269 which you can also know for block-Jacobi, but it works differently and it is only implemented for block 

270 Gauss-Seidel so far. 

271 

272 There is an option to reduce restarts if continued iterations could yield convergence in fewer iterations than 

273 restarting based on an estimate of the contraction factor. 

274 Since often only one or two more iterations suffice, this can boost efficiency of adaptivity significantly. 

275 Notice that the computed step size is not effected. 

276 Be aware that this does not work when Hot Rod is enabled, since that requires us to know the order of the scheme in 

277 more detail. Since we reset to the second to last sweep before moving on, we cannot continue to iterate. 

278 Set the reduced restart up by setting a boolean value for "avoid_restarts" in the parameters for the convergence 

279 controller. 

280 The behaviour in multi-step SDC is not well studied and it is unclear if anything useful happens there. 

281 """ 

282 

283 def setup(self, controller, params, description, **kwargs): 

284 """ 

285 Define default parameters here. 

286 

287 Default parameters are: 

288 - control_order (int): The order relative to other convergence controllers 

289 - beta (float): The safety factor 

290 

291 Args: 

292 controller (pySDC.Controller): The controller 

293 params (dict): The params passed for this specific convergence controller 

294 description (dict): The description object used to instantiate the controller 

295 

296 Returns: 

297 (dict): The updated params dictionary 

298 """ 

299 defaults = { 

300 "embedded_error_flavor": 'standard', 

301 } 

302 return {**defaults, **super().setup(controller, params, description, **kwargs)} 

303 

304 def dependencies(self, controller, description, **kwargs): 

305 """ 

306 Load the embedded error estimator. 

307 

308 Args: 

309 controller (pySDC.Controller): The controller 

310 description (dict): The description object used to instantiate the controller 

311 

312 Returns: 

313 None 

314 """ 

315 from pySDC.implementations.convergence_controller_classes.estimate_embedded_error import EstimateEmbeddedError 

316 

317 super().dependencies(controller, description, **kwargs) 

318 

319 controller.add_convergence_controller( 

320 EstimateEmbeddedError.get_implementation(self.params.embedded_error_flavor, self.params.useMPI), 

321 description=description, 

322 ) 

323 

324 # load contraction factor estimator if necessary 

325 if self.params.get('avoid_restarts'): 

326 from pySDC.implementations.convergence_controller_classes.estimate_contraction_factor import ( 

327 EstimateContractionFactor, 

328 ) 

329 

330 params = {'e_tol': self.params.e_tol} 

331 controller.add_convergence_controller(EstimateContractionFactor, description=description, params=params) 

332 return None 

333 

334 def check_parameters(self, controller, params, description, **kwargs): 

335 """ 

336 Check whether parameters are compatible with whatever assumptions went into the step size functions etc. 

337 For adaptivity, we need to know the order of the scheme. 

338 

339 Args: 

340 controller (pySDC.Controller): The controller 

341 params (dict): The params passed for this specific convergence controller 

342 description (dict): The description object used to instantiate the controller 

343 

344 Returns: 

345 bool: Whether the parameters are compatible 

346 str: The error message 

347 """ 

348 if description["level_params"].get("restol", -1.0) >= 0: 

349 return ( 

350 False, 

351 "Adaptivity needs constant order in time and hence restol in the step parameters has to be \ 

352smaller than 0!", 

353 ) 

354 

355 if controller.params.mssdc_jac: 

356 return ( 

357 False, 

358 "Adaptivity needs the same order on all steps, please activate Gauss-Seidel multistep mode!", 

359 ) 

360 

361 if "e_tol" not in params.keys(): 

362 return ( 

363 False, 

364 "Adaptivity needs a local tolerance! Please pass `e_tol` to the parameters for this convergence controller!", 

365 ) 

366 

367 return True, "" 

368 

369 def get_new_step_size(self, controller, S, **kwargs): 

370 """ 

371 Determine a step size for the next step from an embedded estimate of the local error of the current step. 

372 

373 Args: 

374 controller (pySDC.Controller): The controller 

375 S (pySDC.Step): The current step 

376 

377 Returns: 

378 None 

379 """ 

380 # check if we performed the desired amount of sweeps 

381 if S.status.iter == S.params.maxiter: 

382 L = S.levels[0] 

383 

384 # compute next step size 

385 order = S.status.iter # embedded error estimate is same order as time marching 

386 

387 e_est = self.get_local_error_estimate(controller, S) 

388 L.status.dt_new = self.compute_optimal_step_size( 

389 self.params.beta, L.params.dt, self.params.e_tol, e_est, order 

390 ) 

391 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S) 

392 

393 return None 

394 

395 def get_local_error_estimate(self, controller, S, **kwargs): 

396 """ 

397 Get the embedded error estimate of the finest level of the step. 

398 

399 Args: 

400 controller (pySDC.Controller): The controller 

401 S (pySDC.Step): The current step 

402 

403 Returns: 

404 float: Embedded error estimate 

405 """ 

406 return S.levels[0].status.error_embedded_estimate 

407 

408 

409class AdaptivityRK(Adaptivity): 

410 """ 

411 Adaptivity for Runge-Kutta methods. Basically, we need to change the order in the step size update 

412 """ 

413 

414 def setup(self, controller, params, description, **kwargs): 

415 defaults = {} 

416 defaults['update_order'] = params.get('update_order', description['sweeper_class'].get_update_order()) 

417 return {**defaults, **super().setup(controller, params, description, **kwargs)} 

418 

419 def get_new_step_size(self, controller, S, **kwargs): 

420 """ 

421 Determine a step size for the next step from an embedded estimate of the local error of the current step. 

422 

423 Args: 

424 controller (pySDC.Controller): The controller 

425 S (pySDC.Step): The current step 

426 

427 Returns: 

428 None 

429 """ 

430 # check if we performed the desired amount of sweeps 

431 if S.status.iter == S.params.maxiter: 

432 L = S.levels[0] 

433 

434 # compute next step size 

435 order = self.params.update_order 

436 e_est = self.get_local_error_estimate(controller, S) 

437 L.status.dt_new = self.compute_optimal_step_size( 

438 self.params.beta, L.params.dt, self.params.e_tol, e_est, order 

439 ) 

440 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S) 

441 

442 return None 

443 

444 

445class AdaptivityResidual(AdaptivityBase): 

446 """ 

447 Do adaptivity based on residual. 

448 

449 Since we don't know a correlation between the residual and the error (for nonlinear problems), we employ a simpler 

450 rule to update the step size. Instead of giving a local tolerance that we try to hit as closely as possible, we set 

451 two thresholds for the residual. When we exceed the upper one, we reduce the step size by a factor of 2 and if the 

452 residual falls below the lower threshold, we double the step size. 

453 Please setup these parameters as "e_tol" and "e_tol_low". 

454 """ 

455 

456 def setup(self, controller, params, description, **kwargs): 

457 """ 

458 Define default parameters here. 

459 

460 Default parameters are: 

461 - control_order (int): The order relative to other convergence controllers 

462 - e_tol_low (float): Lower absolute threshold for the residual 

463 - e_tol (float): Upper absolute threshold for the residual 

464 - use_restol (bool): Restart if the residual tolerance was not reached 

465 - max_restarts: Override maximum number of restarts 

466 

467 Args: 

468 controller (pySDC.Controller): The controller 

469 params (dict): The params passed for this specific convergence controller 

470 description (dict): The description object used to instantiate the controller 

471 

472 Returns: 

473 (dict): The updated params dictionary 

474 """ 

475 defaults = { 

476 "control_order": -45, 

477 "e_tol_low": 0, 

478 "e_tol": np.inf, 

479 "use_restol": False, 

480 "max_restarts": 99 if "e_tol_low" in params else None, 

481 "allowed_modifications": ['increase', 'decrease'], # what we are allowed to do with the step size 

482 } 

483 return {**defaults, **params} 

484 

485 def setup_status_variables(self, controller, **kwargs): 

486 """ 

487 Change maximum number of allowed restarts here. 

488 

489 Args: 

490 controller (pySDC.Controller): The controller 

491 

492 Reutrns: 

493 None 

494 """ 

495 from pySDC.implementations.convergence_controller_classes.basic_restarting import BasicRestarting 

496 

497 if self.params.max_restarts is not None: 

498 conv_controllers = controller.convergence_controllers 

499 restart_cont = [me for me in conv_controllers if BasicRestarting in type(me).__bases__] 

500 

501 if len(restart_cont) == 0: 

502 raise NotImplementedError("Please implement override of maximum number of restarts!") 

503 

504 restart_cont[0].params.max_restarts = self.params.max_restarts 

505 return None 

506 

507 def check_parameters(self, controller, params, description, **kwargs): 

508 """ 

509 Check whether parameters are compatible with whatever assumptions went into the step size functions etc. 

510 

511 Args: 

512 controller (pySDC.Controller): The controller 

513 params (dict): The params passed for this specific convergence controller 

514 description (dict): The description object used to instantiate the controller 

515 

516 Returns: 

517 bool: Whether the parameters are compatible 

518 str: The error message 

519 """ 

520 if controller.params.mssdc_jac: 

521 return ( 

522 False, 

523 "Adaptivity needs the same order on all steps, please activate Gauss-Seidel multistep mode!", 

524 ) 

525 

526 return True, "" 

527 

528 def get_new_step_size(self, controller, S, **kwargs): 

529 """ 

530 Determine a step size for the next step. 

531 If we exceed the absolute tolerance of the residual in either direction, we either double or halve the step 

532 size. 

533 

534 Args: 

535 controller (pySDC.Controller): The controller 

536 S (pySDC.Step): The current step 

537 

538 Returns: 

539 None 

540 """ 

541 # check if we performed the desired amount of sweeps 

542 if S.status.iter == S.params.maxiter: 

543 L = S.levels[0] 

544 

545 res = self.get_local_error_estimate(controller, S) 

546 

547 dt_planned = L.status.dt_new if L.status.dt_new is not None else L.params.dt 

548 

549 if ( 

550 res > self.params.e_tol or (res > L.params.restol and self.params.use_restol) 

551 ) and 'decrease' in self.params.allowed_modifications: 

552 L.status.dt_new = min([dt_planned, L.params.dt / 2.0]) 

553 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S) 

554 elif res < self.params.e_tol_low and 'increase' in self.params.allowed_modifications: 

555 L.status.dt_new = max([dt_planned, L.params.dt * 2.0]) 

556 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S) 

557 

558 return None 

559 

560 def get_local_error_estimate(self, controller, S, **kwargs): 

561 """ 

562 Get the residual of the finest level of the step. 

563 

564 Args: 

565 controller (pySDC.Controller): The controller 

566 S (pySDC.Step): The current step 

567 

568 Returns: 

569 float: Embedded error estimate 

570 """ 

571 return S.levels[0].status.residual 

572 

573 

574class AdaptivityCollocation(AdaptivityForConvergedCollocationProblems): 

575 """ 

576 Control the step size via a collocation based estimate of the local error. 

577 The error estimate works by subtracting two solutions to collocation problems with different order. You can 

578 interpolate between collocation methods as much as you want but the adaptive step size selection will always be 

579 based on the last switch of quadrature. 

580 """ 

581 

582 def setup(self, controller, params, description, **kwargs): 

583 """ 

584 Add a default value for control order to the parameters. 

585 

586 Args: 

587 controller (pySDC.Controller): The controller 

588 params (dict): Parameters for the convergence controller 

589 description (dict): The description object used to instantiate the controller 

590 

591 Returns: 

592 dict: Updated parameters 

593 """ 

594 defaults = { 

595 "adaptive_coll_params": {}, 

596 "num_colls": 0, 

597 **super().setup(controller, params, description, **kwargs), 

598 "control_order": 220, 

599 } 

600 

601 for key in defaults['adaptive_coll_params'].keys(): 

602 if type(defaults['adaptive_coll_params'][key]) == list: 

603 defaults['num_colls'] = max([defaults['num_colls'], len(defaults['adaptive_coll_params'][key])]) 

604 

605 if defaults['restart_at_maxiter']: 

606 defaults['maxiter'] = description['step_params'].get('maxiter', 99) * defaults['num_colls'] 

607 

608 return defaults 

609 

610 def setup_status_variables(self, controller, **kwargs): 

611 self.status = Status(['error', 'order']) 

612 self.status.error = [] 

613 self.status.order = [] 

614 

615 def reset_status_variables(self, controller, **kwargs): 

616 self.setup_status_variables(controller, **kwargs) 

617 

618 def dependencies(self, controller, description, **kwargs): 

619 """ 

620 Load the `EstimateEmbeddedErrorCollocation` convergence controller to estimate the local error by switching 

621 between collocation problems between iterations. 

622 

623 Args: 

624 controller (pySDC.Controller): The controller 

625 description (dict): The description object used to instantiate the controller 

626 """ 

627 from pySDC.implementations.convergence_controller_classes.estimate_embedded_error import ( 

628 EstimateEmbeddedErrorCollocation, 

629 ) 

630 

631 super().dependencies(controller, description, **kwargs) 

632 

633 params = {'adaptive_coll_params': self.params.adaptive_coll_params} 

634 controller.add_convergence_controller( 

635 EstimateEmbeddedErrorCollocation, 

636 params=params, 

637 description=description, 

638 ) 

639 

640 def get_convergence(self, controller, S, **kwargs): 

641 return len(self.status.order) == self.params.num_colls 

642 

643 def get_local_error_estimate(self, controller, S, **kwargs): 

644 """ 

645 Get the collocation based embedded error estimate. 

646 

647 Args: 

648 controller (pySDC.Controller): The controller 

649 S (pySDC.Step): The current step 

650 

651 Returns: 

652 float: Embedded error estimate 

653 """ 

654 if len(self.status.error) > 1: 

655 return self.status.error[-1][1] 

656 else: 

657 return 0.0 

658 

659 def post_iteration_processing(self, controller, step, **kwargs): 

660 """ 

661 Get the error estimate and its order if available. 

662 

663 Args: 

664 controller (pySDC.Controller.controller): The controller 

665 step (pySDC.Step.step): The current step 

666 """ 

667 if step.status.done: 

668 lvl = step.levels[0] 

669 self.status.error += [lvl.status.error_embedded_estimate_collocation] 

670 self.status.order += [lvl.sweep.coll.order] 

671 

672 def get_new_step_size(self, controller, S, **kwargs): 

673 if len(self.status.order) == self.params.num_colls: 

674 lvl = S.levels[0] 

675 

676 # compute next step size 

677 order = ( 

678 min(self.status.order[-2::]) + 1 

679 ) # local order of less accurate of the last two collocation problems 

680 e_est = self.get_local_error_estimate(controller, S) 

681 

682 lvl.status.dt_new = self.compute_optimal_step_size( 

683 self.params.beta, lvl.params.dt, self.params.e_tol, e_est, order 

684 ) 

685 self.log(f'Adjusting step size from {lvl.params.dt:.2e} to {lvl.status.dt_new:.2e}', S) 

686 

687 def determine_restart(self, controller, S, **kwargs): 

688 """ 

689 Check if the step wants to be restarted by comparing the estimate of the local error to a preset tolerance 

690 

691 Args: 

692 controller (pySDC.Controller): The controller 

693 S (pySDC.Step): The current step 

694 

695 Returns: 

696 None 

697 """ 

698 if len(self.status.order) == self.params.num_colls: 

699 e_est = self.get_local_error_estimate(controller, S) 

700 if e_est >= self.params.e_tol: 

701 S.status.restart = True 

702 self.log(f"Restarting: e={e_est:.2e} >= e_tol={self.params.e_tol:.2e}", S) 

703 

704 def check_parameters(self, controller, params, description, **kwargs): 

705 """ 

706 Check whether parameters are compatible with whatever assumptions went into the step size functions etc. 

707 For adaptivity, we need to know the order of the scheme. 

708 

709 Args: 

710 controller (pySDC.Controller): The controller 

711 params (dict): The params passed for this specific convergence controller 

712 description (dict): The description object used to instantiate the controller 

713 

714 Returns: 

715 bool: Whether the parameters are compatible 

716 str: The error message 

717 """ 

718 if "e_tol" not in params.keys(): 

719 return ( 

720 False, 

721 "Adaptivity needs a local tolerance! Please pass `e_tol` to the parameters for this convergence controller!", 

722 ) 

723 

724 return True, "" 

725 

726 

727class AdaptivityExtrapolationWithinQ(AdaptivityForConvergedCollocationProblems): 

728 """ 

729 Class to compute time step size adaptively based on error estimate obtained from extrapolation within the quadrature 

730 nodes. 

731 

732 This error estimate depends on solving the collocation problem exactly, so make sure you set a sufficient stopping criterion. 

733 """ 

734 

735 def setup(self, controller, params, description, **kwargs): 

736 from pySDC.implementations.convergence_controller_classes.check_convergence import CheckConvergence 

737 

738 defaults = { 

739 'high_Taylor_order': False, 

740 **params, 

741 } 

742 

743 self.check_convergence = CheckConvergence.check_convergence 

744 return {**defaults, **super().setup(controller, params, description, **kwargs)} 

745 

746 def get_convergence(self, controller, S, **kwargs): 

747 return self.check_convergence(S) 

748 

749 def dependencies(self, controller, description, **kwargs): 

750 """ 

751 Load the error estimator. 

752 

753 Args: 

754 controller (pySDC.Controller): The controller 

755 description (dict): The description object used to instantiate the controller 

756 

757 Returns: 

758 None 

759 """ 

760 from pySDC.implementations.convergence_controller_classes.estimate_extrapolation_error import ( 

761 EstimateExtrapolationErrorWithinQ, 

762 ) 

763 

764 super().dependencies(controller, description, **kwargs) 

765 

766 controller.add_convergence_controller( 

767 EstimateExtrapolationErrorWithinQ, 

768 description=description, 

769 params={'high_Taylor_order': self.params.high_Taylor_order}, 

770 ) 

771 return None 

772 

773 def get_new_step_size(self, controller, S, **kwargs): 

774 """ 

775 Determine a step size for the next step from the error estimate. 

776 

777 Args: 

778 controller (pySDC.Controller): The controller 

779 S (pySDC.Step): The current step 

780 

781 Returns: 

782 None 

783 """ 

784 if self.get_convergence(controller, S, **kwargs): 

785 L = S.levels[0] 

786 

787 # compute next step size 

788 order = L.sweep.coll.num_nodes + 1 if self.params.high_Taylor_order else L.sweep.coll.num_nodes 

789 

790 e_est = self.get_local_error_estimate(controller, S) 

791 L.status.dt_new = self.compute_optimal_step_size( 

792 self.params.beta, L.params.dt, self.params.e_tol, e_est, order 

793 ) 

794 

795 self.log( 

796 f'Error target: {self.params.e_tol:.2e}, error estimate: {e_est:.2e}, update_order: {order}', 

797 S, 

798 level=10, 

799 ) 

800 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S) 

801 

802 return None 

803 

804 def get_local_error_estimate(self, controller, S, **kwargs): 

805 """ 

806 Get the embedded error estimate of the finest level of the step. 

807 

808 Args: 

809 controller (pySDC.Controller): The controller 

810 S (pySDC.Step): The current step 

811 

812 Returns: 

813 float: Embedded error estimate 

814 """ 

815 return S.levels[0].status.error_extrapolation_estimate 

816 

817 

818class AdaptivityPolynomialError(AdaptivityForConvergedCollocationProblems): 

819 """ 

820 Class to compute time step size adaptively based on error estimate obtained from interpolation within the quadrature 

821 nodes. 

822 

823 This error estimate depends on solving the collocation problem exactly, so make sure you set a sufficient stopping criterion. 

824 """ 

825 

826 def setup(self, controller, params, description, **kwargs): 

827 from pySDC.implementations.convergence_controller_classes.check_convergence import CheckConvergence 

828 

829 defaults = { 

830 'control_order': -50, 

831 **super().setup(controller, params, description, **kwargs), 

832 **params, 

833 } 

834 

835 self.check_convergence = CheckConvergence.check_convergence 

836 return defaults 

837 

838 def get_convergence(self, controller, S, **kwargs): 

839 return self.check_convergence(S) 

840 

841 def dependencies(self, controller, description, **kwargs): 

842 """ 

843 Load the error estimator. 

844 

845 Args: 

846 controller (pySDC.Controller): The controller 

847 description (dict): The description object used to instantiate the controller 

848 

849 Returns: 

850 None 

851 """ 

852 from pySDC.implementations.convergence_controller_classes.estimate_polynomial_error import ( 

853 EstimatePolynomialError, 

854 ) 

855 

856 super().dependencies(controller, description, **kwargs) 

857 

858 controller.add_convergence_controller( 

859 EstimatePolynomialError, 

860 description=description, 

861 params={}, 

862 ) 

863 return None 

864 

865 def get_new_step_size(self, controller, S, **kwargs): 

866 """ 

867 Determine a step size for the next step from the error estimate. 

868 

869 Args: 

870 controller (pySDC.Controller): The controller 

871 S (pySDC.Step): The current step 

872 

873 Returns: 

874 None 

875 """ 

876 if self.get_convergence(controller, S, **kwargs): 

877 L = S.levels[0] 

878 

879 # compute next step size 

880 order = L.status.order_embedded_estimate 

881 

882 e_est = self.get_local_error_estimate(controller, S) 

883 L.status.dt_new = self.compute_optimal_step_size( 

884 self.params.beta, L.params.dt, self.params.e_tol, e_est, order 

885 ) 

886 

887 self.log( 

888 f'Error target: {self.params.e_tol:.2e}, error estimate: {e_est:.2e}, update_order: {order}', 

889 S, 

890 level=10, 

891 ) 

892 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S) 

893 

894 return None 

895 

896 def get_local_error_estimate(self, controller, S, **kwargs): 

897 """ 

898 Get the embedded error estimate of the finest level of the step. 

899 

900 Args: 

901 controller (pySDC.Controller): The controller 

902 S (pySDC.Step): The current step 

903 

904 Returns: 

905 float: Embedded error estimate 

906 """ 

907 return S.levels[0].status.error_embedded_estimate