Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
slepc_spectral_transformation.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2018 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 #include <deal.II/lac/slepc_spectral_transformation.h>
17 
18 #ifdef DEAL_II_WITH_SLEPC
19 
20 # include <deal.II/lac/petsc_matrix_base.h>
21 # include <deal.II/lac/slepc_solver.h>
22 
23 # include <petscversion.h>
24 
25 # include <cmath>
26 # include <vector>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace SLEPcWrappers
31 {
32  TransformationBase::TransformationBase(const MPI_Comm &mpi_communicator)
33  {
34  const PetscErrorCode ierr = STCreate(mpi_communicator, &st);
35  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
36  }
37 
39  {
40  if (st != nullptr)
41  {
42  const PetscErrorCode ierr = STDestroy(&st);
43  (void)ierr;
44  AssertNothrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
45  }
46  }
47 
48  void
50  {
51  const PetscErrorCode ierr = STSetMatMode(st, mode);
52  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
53  }
54 
55  void
57  {
58  PetscErrorCode ierr = STSetKSP(st, solver.solver_data->ksp);
59  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
60  }
61 
62  /* ------------------- TransformationShift --------------------- */
63 
65  const double shift_parameter)
66  : shift_parameter(shift_parameter)
67  {}
68 
69  TransformationShift::TransformationShift(const MPI_Comm &mpi_communicator,
70  const AdditionalData &data)
71  : TransformationBase(mpi_communicator)
72  , additional_data(data)
73  {
74  PetscErrorCode ierr = STSetType(st, const_cast<char *>(STSHIFT));
75  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
76 
77  ierr = STSetShift(st, additional_data.shift_parameter);
78  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
79  }
80 
81  /* ---------------- TransformationShiftInvert ------------------ */
82 
84  const double shift_parameter)
85  : shift_parameter(shift_parameter)
86  {}
87 
89  const MPI_Comm & mpi_communicator,
90  const AdditionalData &data)
91  : TransformationBase(mpi_communicator)
92  , additional_data(data)
93  {
94  PetscErrorCode ierr = STSetType(st, const_cast<char *>(STSINVERT));
95  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
96 
97  ierr = STSetShift(st, additional_data.shift_parameter);
98  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
99  }
100 
101  /* --------------- TransformationSpectrumFolding ----------------- */
102 
104  const double shift_parameter)
105  : shift_parameter(shift_parameter)
106  {}
107 
109  const MPI_Comm & mpi_communicator,
110  const AdditionalData &data)
111  : TransformationBase(mpi_communicator)
112  , additional_data(data)
113  {
114 # if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
115  PetscErrorCode ierr = STSetType(st, const_cast<char *>(STFOLD));
116  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
117 
118  ierr = STSetShift(st, additional_data.shift_parameter);
119  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
120 # else
121  // PETSc/SLEPc version must be < 3.5.0.
122  (void)st;
123  Assert((false),
124  ExcMessage(
125  "Folding transformation has been removed in SLEPc 3.5.0 and newer."
126  "You cannot use this transformation anymore."));
127 # endif
128  }
129 
130  /* ------------------- TransformationCayley --------------------- */
131 
133  const double shift_parameter,
134  const double antishift_parameter)
135  : shift_parameter(shift_parameter)
136  , antishift_parameter(antishift_parameter)
137  {}
138 
139  TransformationCayley::TransformationCayley(const MPI_Comm &mpi_communicator,
140  const AdditionalData &data)
141  : TransformationBase(mpi_communicator)
142  , additional_data(data)
143  {
144  PetscErrorCode ierr = STSetType(st, const_cast<char *>(STCAYLEY));
145  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
146 
147  ierr = STSetShift(st, additional_data.shift_parameter);
148  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
149 
150  ierr = STCayleySetAntishift(st, additional_data.antishift_parameter);
151  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
152  }
153 
154 } // namespace SLEPcWrappers
155 
156 DEAL_II_NAMESPACE_CLOSE
157 
158 #endif // DEAL_II_WITH_SLEPC
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1471
void set_solver(const PETScWrappers::SolverBase &solver)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1519
TransformationShiftInvert(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcMessage(std::string arg1)
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:247
#define Assert(cond, exc)
Definition: exceptions.h:1407
TransformationCayley(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
TransformationBase(const MPI_Comm &mpi_communicator)
TransformationShift(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
AdditionalData(const double shift_parameter=0, const double antishift_parameter=0)
static ::ExceptionBase & ExcSLEPcError(int arg1)
TransformationSpectrumFolding(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())