Reference documentation for deal.II version 8.5.1
slepc_spectral_transformation.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2016 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 # include <deal.II/lac/petsc_vector_base.h>
23 # include <deal.II/lac/petsc_vector.h>
24 
25 # include <cmath>
26 # include <vector>
27 
28 # include <petscversion.h>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace SLEPcWrappers
33 {
34  TransformationBase::TransformationBase (const MPI_Comm &mpi_communicator)
35  {
36  const PetscErrorCode ierr = STCreate(mpi_communicator, &st);
37  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
38  }
39 
41  {
42  if (st!=NULL)
43  {
44  const PetscErrorCode ierr = STDestroy(&st);
45  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
46  }
47  }
48 
49  void TransformationBase::set_matrix_mode(const STMatMode mode)
50  {
51  const PetscErrorCode ierr = STSetMatMode(st,mode);
52  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
53  }
54 
56  {
57  PetscErrorCode ierr = STSetKSP(st,solver.solver_data->ksp);
58  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
59  }
60 
61  /* ------------------- TransformationShift --------------------- */
62 
64  AdditionalData (const double shift_parameter)
65  :
66  shift_parameter (shift_parameter)
67  {}
68 
69  TransformationShift::TransformationShift (const MPI_Comm &mpi_communicator,
70  const AdditionalData &data)
71  :
72  TransformationBase(mpi_communicator),
73  additional_data (data)
74  {
75  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STSHIFT));
76  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
77 
78  ierr = STSetShift (st, additional_data.shift_parameter);
79  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
80  }
81 
82  /* ---------------- TransformationShiftInvert ------------------ */
83 
85  AdditionalData (const double shift_parameter)
86  :
87  shift_parameter (shift_parameter)
88  {}
89 
90  TransformationShiftInvert::TransformationShiftInvert (const MPI_Comm &mpi_communicator,
91  const AdditionalData &data)
92  :
93  TransformationBase(mpi_communicator),
94  additional_data (data)
95  {
96 #if DEAL_II_PETSC_VERSION_LT(3,1,0)
97  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STSINV));
98 #else
99  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STSINVERT));
100 #endif
101  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
102 
103  ierr = STSetShift (st, additional_data.shift_parameter);
104  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
105  }
106 
107  /* --------------- TransformationSpectrumFolding ----------------- */
108 
110  AdditionalData (const double shift_parameter)
111  :
112  shift_parameter (shift_parameter)
113  {}
114 
116  const AdditionalData &data)
117  :
118  TransformationBase(mpi_communicator),
119  additional_data (data)
120  {
121 #if DEAL_II_PETSC_VERSION_LT(3,5,0)
122  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STFOLD));
123  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
124 
125  ierr = STSetShift (st, additional_data.shift_parameter);
126  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
127 #else
128  // PETSc/SLEPc version must be < 3.5.0.
129  (void)st;
130  Assert ((false),
131  ExcMessage ("Folding transformation has been removed in SLEPc 3.5.0 and newer."
132  "You cannot use this transformation anymore."));
133 #endif
134  }
135 
136  /* ------------------- TransformationCayley --------------------- */
137 
139  AdditionalData (const double shift_parameter,
140  const double antishift_parameter)
141  :
142  shift_parameter (shift_parameter),
143  antishift_parameter (antishift_parameter)
144  {
145  }
146 
147  TransformationCayley::TransformationCayley (const MPI_Comm &mpi_communicator,
148  const AdditionalData &data)
149  :
150  TransformationBase(mpi_communicator),
151  additional_data (data)
152  {
153  PetscErrorCode ierr = STSetType (st, const_cast<char *>(STCAYLEY));
154  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
155 
156  ierr = STSetShift (st, additional_data.shift_parameter);
157  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
158 
159  ierr = STCayleySetAntishift (st, additional_data.antishift_parameter);
160  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
161  }
162 
163 }
164 
165 DEAL_II_NAMESPACE_CLOSE
166 
167 #endif // DEAL_II_WITH_SLEPC
void set_solver(const PETScWrappers::SolverBase &solver)
std_cxx11::shared_ptr< SolverData > solver_data
Definition: petsc_solver.h:253
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
TransformationShiftInvert(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
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())