Reference documentation for deal.II version GIT 5983d193e2 2023-05-27 16:50:02+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\}}\)
slepc_spectral_transformation.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2019 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 
17 
18 #ifdef DEAL_II_WITH_SLEPC
19 
22 
23 # include <petscversion.h>
24 
25 # include <cmath>
26 # include <vector>
27 
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);
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  // This feature is only in PETSc/SLEPc versions 3.4 or older, which we no
115  // longer support.
116  Assert(false,
117  ExcMessage(
118  "Folding transformation has been removed in SLEPc 3.5.0 and newer."
119  " You cannot use this transformation anymore."));
120  }
121 
122  /* ------------------- TransformationCayley --------------------- */
123 
125  const double shift_parameter,
126  const double antishift_parameter)
127  : shift_parameter(shift_parameter)
128  , antishift_parameter(antishift_parameter)
129  {}
130 
131  TransformationCayley::TransformationCayley(const MPI_Comm mpi_communicator,
132  const AdditionalData &data)
133  : TransformationBase(mpi_communicator)
134  , additional_data(data)
135  {
136  PetscErrorCode ierr = STSetType(st, const_cast<char *>(STCAYLEY));
137  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
138 
139  ierr = STSetShift(st, additional_data.shift_parameter);
140  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
141 
142  ierr = STCayleySetAntishift(st, additional_data.antishift_parameter);
143  AssertThrow(ierr == 0, SolverBase::ExcSLEPcError(ierr));
144  }
145 
146 } // namespace SLEPcWrappers
147 
149 
150 #endif // DEAL_II_WITH_SLEPC
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:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
static ::ExceptionBase & ExcSLEPcError(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1614
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1656
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1703
AdditionalData(const double shift_parameter=0, const double antishift_parameter=0)