Reference documentation for deal.II version 9.3.3
\(\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\}}\)
slepc_spectral_transformation.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2009 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
17#ifndef dealii_slepc_spectral_transformation_h
18# define dealii_slepc_spectral_transformation_h
19
20
21# include <deal.II/base/config.h>
22
23# ifdef DEAL_II_WITH_SLEPC
24
27
28# include <petscksp.h>
29
30# include <slepceps.h>
31
32# include <memory>
33
35
36// Forward declaration
37# ifndef DOXYGEN
38namespace PETScWrappers
39{
40 // forward declarations
41 class SolverBase;
42} // namespace PETScWrappers
43# endif
44
45namespace SLEPcWrappers
46{
47 // forward declaration
48 class SolverBase;
49
77 {
78 protected:
82 TransformationBase(const MPI_Comm &mpi_communicator);
83
84 public:
88 virtual ~TransformationBase();
89
98 void
99 set_matrix_mode(const STMatMode mode);
100
105 void
107
108 protected:
112 ST st;
113
114 // Make the solver class a friend, since it needs to set spectral
115 // transformation object.
116 friend class SolverBase;
117 };
118
125 {
126 public:
131 {
135 AdditionalData(const double shift_parameter = 0);
136
140 const double shift_parameter;
141 };
142
143
147 TransformationShift(const MPI_Comm & mpi_communicator,
148 const AdditionalData &data = AdditionalData());
149
150
151 protected:
156 };
157
165 {
166 public:
171 {
175 AdditionalData(const double shift_parameter = 0);
176
180 const double shift_parameter;
181 };
182
183
187 TransformationShiftInvert(const MPI_Comm & mpi_communicator,
188 const AdditionalData &data = AdditionalData());
189
190 protected:
195
196 // Make the solver class a friend, since it may need to set target
197 // equal the provided shift value.
198 friend class SolverBase;
199 };
200
209 {
210 public:
215 {
219 AdditionalData(const double shift_parameter = 0);
220
224 const double shift_parameter;
225 };
226
227
232 const MPI_Comm & mpi_communicator,
233 const AdditionalData &data = AdditionalData());
234
235 protected:
240 };
241
248 {
249 public:
254 {
258 AdditionalData(const double shift_parameter = 0,
259 const double antishift_parameter = 0);
260
264 const double shift_parameter;
265
270 };
271
272
276 TransformationCayley(const MPI_Comm & mpi_communicator,
277 const AdditionalData &data = AdditionalData());
278
279 protected:
284 };
285
286} // namespace SLEPcWrappers
287
289
290# endif // DEAL_II_WITH_SLEPC
291
292/*-------------------- slepc_spectral_transformation.h ------------------*/
293
294#endif
295
296/*-------------------- slepc_spectral_transformation.h ------------------*/
TransformationBase(const MPI_Comm &mpi_communicator)
void set_solver(const PETScWrappers::SolverBase &solver)
TransformationCayley(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
TransformationShiftInvert(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
TransformationShift(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
TransformationSpectrumFolding(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
AdditionalData(const double shift_parameter=0, const double antishift_parameter=0)