Coverage for pySDC/implementations/convergence_controller_classes/adaptivity.py: 93%
259 statements
« prev ^ index » next coverage.py v7.6.9, created at 2024-12-20 14:51 +0000
« prev ^ index » next coverage.py v7.6.9, created at 2024-12-20 14:51 +0000
1import numpy as np
2from pySDC.core.convergence_controller import ConvergenceController, Status
3from pySDC.implementations.convergence_controller_classes.step_size_limiter import (
4 StepSizeLimiter,
5)
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 """
14 def setup(self, controller, params, description, **kwargs):
15 """
16 Define default parameters here.
18 Default parameters are:
19 - control_order (int): The order relative to other convergence controllers
20 - beta (float): The safety factor
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
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
36 controller.add_hook(LogStepSize)
38 from pySDC.implementations.convergence_controller_classes.check_convergence import CheckConvergence
40 self.communicate_convergence = CheckConvergence.communicate_convergence
42 return {**defaults, **super().setup(controller, params, description, **kwargs)}
44 def dependencies(self, controller, description, **kwargs):
45 """
46 Load step size limiters here, if they are desired.
48 Args:
49 controller (pySDC.Controller): The controller
50 description (dict): The description object used to instantiate the controller
52 Returns:
53 None
54 """
55 step_limiter_keys = ['dt_min', 'dt_max', 'dt_slope_min', 'dt_slope_max', 'dt_rel_min_slope']
56 available_keys = [me for me in step_limiter_keys if me in self.params.__dict__.keys()]
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)
62 if self.params.useMPI:
63 self.prepare_MPI_logical_operations()
65 return None
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.
71 Args:
72 controller (pySDC.Controller): The controller
73 S (pySDC.Step): The current step
75 Returns:
76 None
77 """
78 raise NotImplementedError("Please implement a rule for updating the step size!")
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.
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
94 Returns:
95 float: The optimal step size
96 """
97 return beta * dt * (e_tol / e_est) ** (1.0 / order)
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.
104 Args:
105 controller (pySDC.Controller): The controller
106 S (pySDC.Step): The current step
108 Returns:
109 float: The error estimate
110 """
111 raise NotImplementedError("Please implement a way to get the local error")
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
117 Args:
118 controller (pySDC.Controller): The controller
119 S (pySDC.Step): The current step
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
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)
155 return None
158class AdaptivityForConvergedCollocationProblems(AdaptivityBase):
159 def dependencies(self, controller, description, **kwargs):
160 """
161 Load interpolation between restarts.
163 Args:
164 controller (pySDC.Controller): The controller
165 description (dict): The description object used to instantiate the controller
167 Returns:
168 None
169 """
170 super().dependencies(controller, description, **kwargs)
172 if self.params.interpolate_between_restarts:
173 from pySDC.implementations.convergence_controller_classes.interpolate_between_restarts import (
174 InterpolateBetweenRestarts,
175 )
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
186 def get_convergence(self, controller, S, **kwargs):
187 raise NotImplementedError("Please implement a way to check if the collocation problem is converged!")
189 def setup(self, controller, params, description, **kwargs):
190 """
191 Add a default value for control order to the parameters.
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
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])
221 if defaults['restart_at_maxiter']:
222 defaults['maxiter'] = description['step_params'].get('maxiter', 99)
224 self.res_last_iter = np.inf
226 return defaults
228 def determine_restart(self, controller, S, **kwargs):
229 if self.get_convergence(controller, S, **kwargs):
230 self.res_last_iter = np.inf
232 L = S.levels[0]
233 e_tol_converged = (
234 L.status.increment < L.params.e_tol if (L.params.get('e_tol') and L.status.get('increment')) else False
235 )
237 if (
238 self.params.restart_at_maxiter
239 and S.levels[0].status.residual > S.levels[0].params.restol
240 and not e_tol_converged
241 ):
242 self.trigger_restart_upon_nonconvergence(S)
243 elif self.get_local_error_estimate(controller, S, **kwargs) > self.params.e_tol:
244 S.status.restart = True
245 elif (
246 S.status.time_size == 1
247 and self.res_last_iter < S.levels[0].status.residual
248 and S.status.iter > 0
249 and self.params.abort_at_growing_residual
250 ):
251 self.trigger_restart_upon_nonconvergence(S)
252 elif S.levels[0].status.residual > self.params.residual_max_tol:
253 self.trigger_restart_upon_nonconvergence(S)
255 if self.params.useMPI:
256 self.communicate_convergence(self, controller, S, **kwargs)
258 self.res_last_iter = S.levels[0].status.residual * 1.0
260 def trigger_restart_upon_nonconvergence(self, S):
261 S.status.restart = True
262 S.status.force_done = True
263 for L in S.levels:
264 L.status.dt_new = L.params.dt / self.params.factor_if_not_converged
265 self.log(
266 f'Collocation problem not converged. Reducing step size to {L.status.dt_new:.2e}',
267 S,
268 )
269 if self.params.interpolate_between_restarts:
270 self.interpolator.status.skip_interpolation = True
273class Adaptivity(AdaptivityBase):
274 """
275 Class to compute time step size adaptively based on embedded error estimate.
277 We have a version working in non-MPI pipelined SDC, but Adaptivity requires you to know the order of the scheme,
278 which you can also know for block-Jacobi, but it works differently and it is only implemented for block
279 Gauss-Seidel so far.
281 There is an option to reduce restarts if continued iterations could yield convergence in fewer iterations than
282 restarting based on an estimate of the contraction factor.
283 Since often only one or two more iterations suffice, this can boost efficiency of adaptivity significantly.
284 Notice that the computed step size is not effected.
285 Be aware that this does not work when Hot Rod is enabled, since that requires us to know the order of the scheme in
286 more detail. Since we reset to the second to last sweep before moving on, we cannot continue to iterate.
287 Set the reduced restart up by setting a boolean value for "avoid_restarts" in the parameters for the convergence
288 controller.
289 The behaviour in multi-step SDC is not well studied and it is unclear if anything useful happens there.
290 """
292 def setup(self, controller, params, description, **kwargs):
293 """
294 Define default parameters here.
296 Default parameters are:
297 - control_order (int): The order relative to other convergence controllers
298 - beta (float): The safety factor
300 Args:
301 controller (pySDC.Controller): The controller
302 params (dict): The params passed for this specific convergence controller
303 description (dict): The description object used to instantiate the controller
305 Returns:
306 (dict): The updated params dictionary
307 """
308 defaults = {
309 "embedded_error_flavor": 'standard',
310 }
311 return {**defaults, **super().setup(controller, params, description, **kwargs)}
313 def dependencies(self, controller, description, **kwargs):
314 """
315 Load the embedded error estimator.
317 Args:
318 controller (pySDC.Controller): The controller
319 description (dict): The description object used to instantiate the controller
321 Returns:
322 None
323 """
324 from pySDC.implementations.convergence_controller_classes.estimate_embedded_error import EstimateEmbeddedError
326 super().dependencies(controller, description, **kwargs)
328 controller.add_convergence_controller(
329 EstimateEmbeddedError.get_implementation(self.params.embedded_error_flavor, self.params.useMPI),
330 description=description,
331 )
333 # load contraction factor estimator if necessary
334 if self.params.get('avoid_restarts'):
335 from pySDC.implementations.convergence_controller_classes.estimate_contraction_factor import (
336 EstimateContractionFactor,
337 )
339 params = {'e_tol': self.params.e_tol}
340 controller.add_convergence_controller(EstimateContractionFactor, description=description, params=params)
341 return None
343 def check_parameters(self, controller, params, description, **kwargs):
344 """
345 Check whether parameters are compatible with whatever assumptions went into the step size functions etc.
346 For adaptivity, we need to know the order of the scheme.
348 Args:
349 controller (pySDC.Controller): The controller
350 params (dict): The params passed for this specific convergence controller
351 description (dict): The description object used to instantiate the controller
353 Returns:
354 bool: Whether the parameters are compatible
355 str: The error message
356 """
357 if description["level_params"].get("restol", -1.0) >= 0:
358 return (
359 False,
360 "Adaptivity needs constant order in time and hence restol in the step parameters has to be \
361smaller than 0!",
362 )
364 if controller.params.mssdc_jac:
365 return (
366 False,
367 "Adaptivity needs the same order on all steps, please activate Gauss-Seidel multistep mode!",
368 )
370 if "e_tol" not in params.keys():
371 return (
372 False,
373 "Adaptivity needs a local tolerance! Please pass `e_tol` to the parameters for this convergence controller!",
374 )
376 return True, ""
378 def get_new_step_size(self, controller, S, **kwargs):
379 """
380 Determine a step size for the next step from an embedded estimate of the local error of the current step.
382 Args:
383 controller (pySDC.Controller): The controller
384 S (pySDC.Step): The current step
386 Returns:
387 None
388 """
389 # check if we performed the desired amount of sweeps
390 if S.status.iter == S.params.maxiter:
391 L = S.levels[0]
393 # compute next step size
394 order = S.status.iter # embedded error estimate is same order as time marching
396 e_est = self.get_local_error_estimate(controller, S)
397 L.status.dt_new = self.compute_optimal_step_size(
398 self.params.beta, L.params.dt, self.params.e_tol, e_est, order
399 )
400 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S)
402 return None
404 def get_local_error_estimate(self, controller, S, **kwargs):
405 """
406 Get the embedded error estimate of the finest level of the step.
408 Args:
409 controller (pySDC.Controller): The controller
410 S (pySDC.Step): The current step
412 Returns:
413 float: Embedded error estimate
414 """
415 return S.levels[0].status.error_embedded_estimate
418class AdaptivityRK(Adaptivity):
419 """
420 Adaptivity for Runge-Kutta methods. Basically, we need to change the order in the step size update
421 """
423 def setup(self, controller, params, description, **kwargs):
424 defaults = {}
425 defaults['update_order'] = params.get('update_order', description['sweeper_class'].get_update_order())
426 return {**defaults, **super().setup(controller, params, description, **kwargs)}
428 def get_new_step_size(self, controller, S, **kwargs):
429 """
430 Determine a step size for the next step from an embedded estimate of the local error of the current step.
432 Args:
433 controller (pySDC.Controller): The controller
434 S (pySDC.Step): The current step
436 Returns:
437 None
438 """
439 # check if we performed the desired amount of sweeps
440 if S.status.iter == S.params.maxiter:
441 L = S.levels[0]
443 # compute next step size
444 order = self.params.update_order
445 e_est = self.get_local_error_estimate(controller, S)
446 L.status.dt_new = self.compute_optimal_step_size(
447 self.params.beta, L.params.dt, self.params.e_tol, e_est, order
448 )
449 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S)
451 return None
454class AdaptivityResidual(AdaptivityBase):
455 """
456 Do adaptivity based on residual.
458 Since we don't know a correlation between the residual and the error (for nonlinear problems), we employ a simpler
459 rule to update the step size. Instead of giving a local tolerance that we try to hit as closely as possible, we set
460 two thresholds for the residual. When we exceed the upper one, we reduce the step size by a factor of 2 and if the
461 residual falls below the lower threshold, we double the step size.
462 Please setup these parameters as "e_tol" and "e_tol_low".
463 """
465 def setup(self, controller, params, description, **kwargs):
466 """
467 Define default parameters here.
469 Default parameters are:
470 - control_order (int): The order relative to other convergence controllers
471 - e_tol_low (float): Lower absolute threshold for the residual
472 - e_tol (float): Upper absolute threshold for the residual
473 - use_restol (bool): Restart if the residual tolerance was not reached
474 - max_restarts: Override maximum number of restarts
476 Args:
477 controller (pySDC.Controller): The controller
478 params (dict): The params passed for this specific convergence controller
479 description (dict): The description object used to instantiate the controller
481 Returns:
482 (dict): The updated params dictionary
483 """
484 defaults = {
485 "control_order": -45,
486 "e_tol_low": 0,
487 "e_tol": np.inf,
488 "use_restol": False,
489 "max_restarts": 99 if "e_tol_low" in params else None,
490 "allowed_modifications": ['increase', 'decrease'], # what we are allowed to do with the step size
491 }
492 return {**defaults, **params}
494 def setup_status_variables(self, controller, **kwargs):
495 """
496 Change maximum number of allowed restarts here.
498 Args:
499 controller (pySDC.Controller): The controller
501 Reutrns:
502 None
503 """
504 from pySDC.implementations.convergence_controller_classes.basic_restarting import BasicRestarting
506 if self.params.max_restarts is not None:
507 conv_controllers = controller.convergence_controllers
508 restart_cont = [me for me in conv_controllers if BasicRestarting in type(me).__bases__]
510 if len(restart_cont) == 0:
511 raise NotImplementedError("Please implement override of maximum number of restarts!")
513 restart_cont[0].params.max_restarts = self.params.max_restarts
514 return None
516 def check_parameters(self, controller, params, description, **kwargs):
517 """
518 Check whether parameters are compatible with whatever assumptions went into the step size functions etc.
520 Args:
521 controller (pySDC.Controller): The controller
522 params (dict): The params passed for this specific convergence controller
523 description (dict): The description object used to instantiate the controller
525 Returns:
526 bool: Whether the parameters are compatible
527 str: The error message
528 """
529 if controller.params.mssdc_jac:
530 return (
531 False,
532 "Adaptivity needs the same order on all steps, please activate Gauss-Seidel multistep mode!",
533 )
535 return True, ""
537 def get_new_step_size(self, controller, S, **kwargs):
538 """
539 Determine a step size for the next step.
540 If we exceed the absolute tolerance of the residual in either direction, we either double or halve the step
541 size.
543 Args:
544 controller (pySDC.Controller): The controller
545 S (pySDC.Step): The current step
547 Returns:
548 None
549 """
550 # check if we performed the desired amount of sweeps
551 if S.status.iter == S.params.maxiter:
552 L = S.levels[0]
554 res = self.get_local_error_estimate(controller, S)
556 dt_planned = L.status.dt_new if L.status.dt_new is not None else L.params.dt
558 if (
559 res > self.params.e_tol or (res > L.params.restol and self.params.use_restol)
560 ) and 'decrease' in self.params.allowed_modifications:
561 L.status.dt_new = min([dt_planned, L.params.dt / 2.0])
562 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S)
563 elif res < self.params.e_tol_low and 'increase' in self.params.allowed_modifications:
564 L.status.dt_new = max([dt_planned, L.params.dt * 2.0])
565 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S)
567 return None
569 def get_local_error_estimate(self, controller, S, **kwargs):
570 """
571 Get the residual of the finest level of the step.
573 Args:
574 controller (pySDC.Controller): The controller
575 S (pySDC.Step): The current step
577 Returns:
578 float: Embedded error estimate
579 """
580 return S.levels[0].status.residual
583class AdaptivityCollocation(AdaptivityForConvergedCollocationProblems):
584 """
585 Control the step size via a collocation based estimate of the local error.
586 The error estimate works by subtracting two solutions to collocation problems with different order. You can
587 interpolate between collocation methods as much as you want but the adaptive step size selection will always be
588 based on the last switch of quadrature.
589 """
591 def setup(self, controller, params, description, **kwargs):
592 """
593 Add a default value for control order to the parameters.
595 Args:
596 controller (pySDC.Controller): The controller
597 params (dict): Parameters for the convergence controller
598 description (dict): The description object used to instantiate the controller
600 Returns:
601 dict: Updated parameters
602 """
603 defaults = {
604 "adaptive_coll_params": {},
605 "num_colls": 0,
606 **super().setup(controller, params, description, **kwargs),
607 "control_order": 220,
608 }
610 for key in defaults['adaptive_coll_params'].keys():
611 if type(defaults['adaptive_coll_params'][key]) == list:
612 defaults['num_colls'] = max([defaults['num_colls'], len(defaults['adaptive_coll_params'][key])])
614 if defaults['restart_at_maxiter']:
615 defaults['maxiter'] = description['step_params'].get('maxiter', 99) * defaults['num_colls']
617 return defaults
619 def setup_status_variables(self, controller, **kwargs):
620 self.status = Status(['error', 'order'])
621 self.status.error = []
622 self.status.order = []
624 def reset_status_variables(self, controller, **kwargs):
625 self.setup_status_variables(controller, **kwargs)
627 def dependencies(self, controller, description, **kwargs):
628 """
629 Load the `EstimateEmbeddedErrorCollocation` convergence controller to estimate the local error by switching
630 between collocation problems between iterations.
632 Args:
633 controller (pySDC.Controller): The controller
634 description (dict): The description object used to instantiate the controller
635 """
636 from pySDC.implementations.convergence_controller_classes.estimate_embedded_error import (
637 EstimateEmbeddedErrorCollocation,
638 )
640 super().dependencies(controller, description, **kwargs)
642 params = {'adaptive_coll_params': self.params.adaptive_coll_params}
643 controller.add_convergence_controller(
644 EstimateEmbeddedErrorCollocation,
645 params=params,
646 description=description,
647 )
649 def get_convergence(self, controller, S, **kwargs):
650 return len(self.status.order) == self.params.num_colls
652 def get_local_error_estimate(self, controller, S, **kwargs):
653 """
654 Get the collocation based embedded error estimate.
656 Args:
657 controller (pySDC.Controller): The controller
658 S (pySDC.Step): The current step
660 Returns:
661 float: Embedded error estimate
662 """
663 if len(self.status.error) > 1:
664 return self.status.error[-1][1]
665 else:
666 return 0.0
668 def post_iteration_processing(self, controller, step, **kwargs):
669 """
670 Get the error estimate and its order if available.
672 Args:
673 controller (pySDC.Controller.controller): The controller
674 step (pySDC.Step.step): The current step
675 """
676 if step.status.done:
677 lvl = step.levels[0]
678 self.status.error += [lvl.status.error_embedded_estimate_collocation]
679 self.status.order += [lvl.sweep.coll.order]
681 def get_new_step_size(self, controller, S, **kwargs):
682 if len(self.status.order) == self.params.num_colls:
683 lvl = S.levels[0]
685 # compute next step size
686 order = (
687 min(self.status.order[-2::]) + 1
688 ) # local order of less accurate of the last two collocation problems
689 e_est = self.get_local_error_estimate(controller, S)
691 lvl.status.dt_new = self.compute_optimal_step_size(
692 self.params.beta, lvl.params.dt, self.params.e_tol, e_est, order
693 )
694 self.log(f'Adjusting step size from {lvl.params.dt:.2e} to {lvl.status.dt_new:.2e}', S)
696 def determine_restart(self, controller, S, **kwargs):
697 """
698 Check if the step wants to be restarted by comparing the estimate of the local error to a preset tolerance
700 Args:
701 controller (pySDC.Controller): The controller
702 S (pySDC.Step): The current step
704 Returns:
705 None
706 """
707 if len(self.status.order) == self.params.num_colls:
708 e_est = self.get_local_error_estimate(controller, S)
709 if e_est >= self.params.e_tol:
710 S.status.restart = True
711 self.log(f"Restarting: e={e_est:.2e} >= e_tol={self.params.e_tol:.2e}", S)
713 def check_parameters(self, controller, params, description, **kwargs):
714 """
715 Check whether parameters are compatible with whatever assumptions went into the step size functions etc.
716 For adaptivity, we need to know the order of the scheme.
718 Args:
719 controller (pySDC.Controller): The controller
720 params (dict): The params passed for this specific convergence controller
721 description (dict): The description object used to instantiate the controller
723 Returns:
724 bool: Whether the parameters are compatible
725 str: The error message
726 """
727 if "e_tol" not in params.keys():
728 return (
729 False,
730 "Adaptivity needs a local tolerance! Please pass `e_tol` to the parameters for this convergence controller!",
731 )
733 return True, ""
736class AdaptivityExtrapolationWithinQ(AdaptivityForConvergedCollocationProblems):
737 """
738 Class to compute time step size adaptively based on error estimate obtained from extrapolation within the quadrature
739 nodes.
741 This error estimate depends on solving the collocation problem exactly, so make sure you set a sufficient stopping criterion.
742 """
744 def setup(self, controller, params, description, **kwargs):
745 from pySDC.implementations.convergence_controller_classes.check_convergence import CheckConvergence
747 defaults = {
748 'high_Taylor_order': False,
749 **params,
750 }
752 self.check_convergence = CheckConvergence.check_convergence
753 return {**defaults, **super().setup(controller, params, description, **kwargs)}
755 def get_convergence(self, controller, S, **kwargs):
756 return self.check_convergence(S)
758 def dependencies(self, controller, description, **kwargs):
759 """
760 Load the error estimator.
762 Args:
763 controller (pySDC.Controller): The controller
764 description (dict): The description object used to instantiate the controller
766 Returns:
767 None
768 """
769 from pySDC.implementations.convergence_controller_classes.estimate_extrapolation_error import (
770 EstimateExtrapolationErrorWithinQ,
771 )
773 super().dependencies(controller, description, **kwargs)
775 controller.add_convergence_controller(
776 EstimateExtrapolationErrorWithinQ,
777 description=description,
778 params={'high_Taylor_order': self.params.high_Taylor_order},
779 )
780 return None
782 def get_new_step_size(self, controller, S, **kwargs):
783 """
784 Determine a step size for the next step from the error estimate.
786 Args:
787 controller (pySDC.Controller): The controller
788 S (pySDC.Step): The current step
790 Returns:
791 None
792 """
793 if self.get_convergence(controller, S, **kwargs):
794 L = S.levels[0]
796 # compute next step size
797 order = L.sweep.coll.num_nodes + 1 if self.params.high_Taylor_order else L.sweep.coll.num_nodes
799 e_est = self.get_local_error_estimate(controller, S)
800 L.status.dt_new = self.compute_optimal_step_size(
801 self.params.beta, L.params.dt, self.params.e_tol, e_est, order
802 )
804 self.log(
805 f'Error target: {self.params.e_tol:.2e}, error estimate: {e_est:.2e}, update_order: {order}',
806 S,
807 level=10,
808 )
809 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S)
811 return None
813 def get_local_error_estimate(self, controller, S, **kwargs):
814 """
815 Get the embedded error estimate of the finest level of the step.
817 Args:
818 controller (pySDC.Controller): The controller
819 S (pySDC.Step): The current step
821 Returns:
822 float: Embedded error estimate
823 """
824 return S.levels[0].status.error_extrapolation_estimate
827class AdaptivityPolynomialError(AdaptivityForConvergedCollocationProblems):
828 """
829 Class to compute time step size adaptively based on error estimate obtained from interpolation within the quadrature
830 nodes.
832 This error estimate depends on solving the collocation problem exactly, so make sure you set a sufficient stopping criterion.
833 """
835 def setup(self, controller, params, description, **kwargs):
836 from pySDC.implementations.convergence_controller_classes.check_convergence import CheckConvergence
838 defaults = {
839 'control_order': -50,
840 **super().setup(controller, params, description, **kwargs),
841 **params,
842 }
844 self.check_convergence = CheckConvergence.check_convergence
845 return defaults
847 def get_convergence(self, controller, S, **kwargs):
848 return self.check_convergence(S)
850 def dependencies(self, controller, description, **kwargs):
851 """
852 Load the error estimator.
854 Args:
855 controller (pySDC.Controller): The controller
856 description (dict): The description object used to instantiate the controller
858 Returns:
859 None
860 """
861 from pySDC.implementations.convergence_controller_classes.estimate_polynomial_error import (
862 EstimatePolynomialError,
863 )
865 super().dependencies(controller, description, **kwargs)
867 controller.add_convergence_controller(
868 EstimatePolynomialError,
869 description=description,
870 params={},
871 )
872 return None
874 def get_new_step_size(self, controller, S, **kwargs):
875 """
876 Determine a step size for the next step from the error estimate.
878 Args:
879 controller (pySDC.Controller): The controller
880 S (pySDC.Step): The current step
882 Returns:
883 None
884 """
885 if self.get_convergence(controller, S, **kwargs):
886 L = S.levels[0]
888 # compute next step size
889 order = L.status.order_embedded_estimate
891 e_est = self.get_local_error_estimate(controller, S)
892 L.status.dt_new = self.compute_optimal_step_size(
893 self.params.beta, L.params.dt, self.params.e_tol, e_est, order
894 )
896 self.log(
897 f'Error target: {self.params.e_tol:.2e}, error estimate: {e_est:.2e}, update_order: {order}',
898 S,
899 level=10,
900 )
901 self.log(f'Adjusting step size from {L.params.dt:.2e} to {L.status.dt_new:.2e}', S)
903 return None
905 def get_local_error_estimate(self, controller, S, **kwargs):
906 """
907 Get the embedded error estimate of the finest level of the step.
909 Args:
910 controller (pySDC.Controller): The controller
911 S (pySDC.Step): The current step
913 Returns:
914 float: Embedded error estimate
915 """
916 return S.levels[0].status.error_embedded_estimate