Reference documentation for deal.II version GIT relicensing-1321-g19b0133728 2024-07-26 07:50:01+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
petsc_ts.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2023 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_petsc_ts_h
16#define dealii_petsc_ts_h
17
18#include <deal.II/base/config.h>
19
20#ifdef DEAL_II_WITH_PETSC
21# include <deal.II/base/mpi.h>
24
28
29# include <petscts.h>
30
31# include <exception>
32
33# if defined(DEAL_II_HAVE_CXX20)
34# include <concepts>
35# endif
36
38
39namespace PETScWrappers
40{
45 {
46 public:
52 using real_type = PetscReal;
53
89 // Running parameters
90 const std::string &options_prefix = "",
91 const std::string &ts_type = "",
92 const real_type initial_time = 0.0,
93 const real_type final_time = 0.0,
94 const real_type initial_step_size = 0.0,
95 const int max_steps = -1,
96 const bool match_step = false,
97 const bool restart_if_remesh = false,
98 // Error parameters
99 const std::string &ts_adapt_type = "none",
100 const real_type minimum_step_size = -1.0,
101 const real_type maximum_step_size = -1.0,
102 const real_type absolute_tolerance = -1.0,
103 const real_type relative_tolerance = -1.0,
104 const bool ignore_algebraic_lte = true)
119 {}
120
124 void
126
130 std::string options_prefix;
131
135 std::string ts_type;
136
141
146
153
160
165
170
174 std::string ts_adapt_type;
175
182
189
196
203
208 };
209
326 template <typename VectorType = PETScWrappers::VectorBase,
327 typename PMatrixType = PETScWrappers::MatrixBase,
328 typename AMatrixType = PMatrixType>
329# if defined(DEAL_II_HAVE_CXX20) && \
330 !defined(DEAL_II_DOXYGEN_DO_NOT_PARSE_REQUIRES_CLAUSES)
332 std::constructible_from<VectorType, Vec>) &&
334 std::constructible_from<PMatrixType, Mat>) &&
336 std::constructible_from<AMatrixType, Mat>))
337# endif
339 {
340 public:
346 using real_type = PetscReal;
347
352 const MPI_Comm mpi_comm = PETSC_COMM_WORLD);
353
358
364 operator TS() const;
365
369 TS
371
377
381 void
383
388 void
390
403 void
404 set_matrix(PMatrixType &P);
405
411 void
412 set_matrices(AMatrixType &A, PMatrixType &P);
413
419
425
429 unsigned int
431
439 unsigned int
440 solve(VectorType &y);
441
451 unsigned int
452 solve(VectorType &y, PMatrixType &P);
453
464 unsigned int
465 solve(VectorType &y, AMatrixType &A, PMatrixType &P);
466
475 std::function<void(const real_type t,
476 const VectorType &y,
477 const VectorType &y_dot,
478 VectorType &res)>
480
495 std::function<void(const real_type t,
496 const VectorType &y,
497 const VectorType &y_dot,
498 const real_type alpha,
499 AMatrixType &A,
500 PMatrixType &P)>
502
511 std::function<void(const real_type t, const VectorType &y, VectorType &res)>
513
523 std::function<void(const real_type t,
524 const VectorType &y,
525 AMatrixType &A,
526 PMatrixType &P)>
528
540 std::function<void(const real_type t,
541 const VectorType &y,
542 const unsigned int step_number)>
544
564 std::function<void(const real_type t,
565 const VectorType &y,
566 const VectorType &ydot,
567 const real_type alpha)>
569
581 std::function<void(const VectorType &src, VectorType &dst)>
583
615
621 std::function<void(const real_type t, VectorType &y)> distribute;
622
641 std::function<void(const real_type t, VectorType &y)>
643
652 std::function<void(const real_type t,
653 const unsigned int step,
654 const VectorType &y,
655 bool &resize)>
657
677 std::function<
678 bool(const real_type t, const unsigned int step, const VectorType &y)>
680
686 std::function<void(const std::vector<VectorType> &all_in,
687 std::vector<VectorType> &all_out)>
689
705 std::function<void(const real_type t,
706 const std::vector<VectorType> &all_in,
707 std::vector<VectorType> &all_out)>
709
710 private:
714 TS ts;
715
721
726
732
737
742
748 mutable std::exception_ptr pending_exception;
749
754
762 void
764
772 void
773 setup_algebraic_constraints(const VectorType &y);
774 };
775
776} // namespace PETScWrappers
777
779
780#endif // DEAL_II_WITH_PETSC
781
782#endif
void add_parameters(ParameterHandler &prm)
void reinit(const TimeStepperData &data)
std::exception_ptr pending_exception
Definition petsc_ts.h:748
unsigned int solve(VectorType &y, PMatrixType &P)
std::function< void(const real_type t, const unsigned int step, const VectorType &y, bool &resize)> decide_for_coarsening_and_refinement
Definition petsc_ts.h:656
void setup_algebraic_constraints(const VectorType &y)
std::function< void(const real_type t, const VectorType &y, const VectorType &y_dot, VectorType &res)> implicit_function
Definition petsc_ts.h:479
MPI_Comm get_mpi_communicator() const
std::function< void(const real_type t, VectorType &y)> update_constrained_components
Definition petsc_ts.h:642
std::function< void(const real_type t, const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out)> transfer_solution_vectors_to_new_mesh
Definition petsc_ts.h:708
std::function< void(const real_type t, const VectorType &y, AMatrixType &A, PMatrixType &P)> explicit_jacobian
Definition petsc_ts.h:527
std::function< bool(const real_type t, const unsigned int step, const VectorType &y)> decide_and_prepare_for_remeshing
Definition petsc_ts.h:679
unsigned int get_step_number()
TimeStepper(const TimeStepperData &data=TimeStepperData(), const MPI_Comm mpi_comm=PETSC_COMM_WORLD)
SmartPointer< AMatrixType, TimeStepper > A
Definition petsc_ts.h:719
unsigned int solve(VectorType &y, AMatrixType &A, PMatrixType &P)
std::function< void(const real_type t, const VectorType &y, const unsigned int step_number)> monitor
Definition petsc_ts.h:543
std::function< void(const real_type t, const VectorType &y, const VectorType &y_dot, const real_type alpha, AMatrixType &A, PMatrixType &P)> implicit_jacobian
Definition petsc_ts.h:501
std::function< IndexSet()> algebraic_components
Definition petsc_ts.h:614
std::function< void(const real_type t, VectorType &y)> distribute
Definition petsc_ts.h:621
std::function< void(const real_type t, const VectorType &y, const VectorType &ydot, const real_type alpha)> setup_jacobian
Definition petsc_ts.h:568
SmartPointer< PMatrixType, TimeStepper > P
Definition petsc_ts.h:720
std::function< void(const real_type t, const VectorType &y, VectorType &res)> explicit_function
Definition petsc_ts.h:512
unsigned int solve(VectorType &y)
std::function< void(const VectorType &src, VectorType &dst)> solve_with_jacobian
Definition petsc_ts.h:582
void set_matrices(AMatrixType &A, PMatrixType &P)
PreconditionShell solve_with_jacobian_pc
Definition petsc_ts.h:725
std::function< void(const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out)> interpolate
Definition petsc_ts.h:688
void set_matrix(PMatrixType &P)
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
TimeStepperData(const std::string &options_prefix="", const std::string &ts_type="", const real_type initial_time=0.0, const real_type final_time=0.0, const real_type initial_step_size=0.0, const int max_steps=-1, const bool match_step=false, const bool restart_if_remesh=false, const std::string &ts_adapt_type="none", const real_type minimum_step_size=-1.0, const real_type maximum_step_size=-1.0, const real_type absolute_tolerance=-1.0, const real_type relative_tolerance=-1.0, const bool ignore_algebraic_lte=true)
Definition petsc_ts.h:88