Reference documentation for deal.II version 9.0.0
slepc_spectral_transformation.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2017 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 at
12 // the top level of the deal.II distribution.
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/slepc_solver.h>
21 # include <deal.II/lac/petsc_matrix_base.h>
22 
23 # include <cmath>
24 # include <vector>
25 
26 # include <petscversion.h>
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 TransformationBase::set_matrix_mode(const STMatMode mode)
49  {
50  const PetscErrorCode ierr = STSetMatMode(st,mode);
51  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
52  }
53 
55  {
56  PetscErrorCode ierr = STSetKSP(st,solver.solver_data->ksp);
57  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
58  }
59 
60  /* ------------------- TransformationShift --------------------- */
61 
63  AdditionalData (const double shift_parameter)
64  :
65  shift_parameter (shift_parameter)
66  {}
67 
68  TransformationShift::TransformationShift (const MPI_Comm &mpi_communicator,
69  const AdditionalData &data)
70  :
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  AdditionalData (const double shift_parameter)
85  :
86  shift_parameter (shift_parameter)
87  {}
88 
89  TransformationShiftInvert::TransformationShiftInvert (const MPI_Comm &mpi_communicator,
90  const AdditionalData &data)
91  :
92  TransformationBase(mpi_communicator),
93  additional_data (data)
94  {
95  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STSINVERT));
96  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
97 
98  ierr = STSetShift (st, additional_data.shift_parameter);
99  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
100  }
101 
102  /* --------------- TransformationSpectrumFolding ----------------- */
103 
105  AdditionalData (const double shift_parameter)
106  :
107  shift_parameter (shift_parameter)
108  {}
109 
111  const AdditionalData &data)
112  :
113  TransformationBase(mpi_communicator),
114  additional_data (data)
115  {
116 #if DEAL_II_PETSC_VERSION_LT(3,5,0)
117  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STFOLD));
118  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
119 
120  ierr = STSetShift (st, additional_data.shift_parameter);
121  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
122 #else
123  // PETSc/SLEPc version must be < 3.5.0.
124  (void)st;
125  Assert ((false),
126  ExcMessage ("Folding transformation has been removed in SLEPc 3.5.0 and newer."
127  "You cannot use this transformation anymore."));
128 #endif
129  }
130 
131  /* ------------------- TransformationCayley --------------------- */
132 
134  AdditionalData (const double shift_parameter,
135  const double antishift_parameter)
136  :
137  shift_parameter (shift_parameter),
138  antishift_parameter (antishift_parameter)
139  {
140  }
141 
142  TransformationCayley::TransformationCayley (const MPI_Comm &mpi_communicator,
143  const AdditionalData &data)
144  :
145  TransformationBase(mpi_communicator),
146  additional_data (data)
147  {
148  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STCAYLEY));
149  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
150 
151  ierr = STSetShift (st, additional_data.shift_parameter);
152  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
153 
154  ierr = STCayleySetAntishift (st, additional_data.antishift_parameter);
155  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
156  }
157 
158 }
159 
160 DEAL_II_NAMESPACE_CLOSE
161 
162 #endif // DEAL_II_WITH_SLEPC
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1183
void set_solver(const PETScWrappers::SolverBase &solver)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
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:243
#define Assert(cond, exc)
Definition: exceptions.h:1142
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())