Reference documentation for deal.II version 8.4.2
slepc_spectral_transformation.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2015 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  int ierr = STCreate(mpi_communicator, &st);
37  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
38  }
39 
41  {
42  if (st!=NULL)
43  {
44  int ierr = STDestroy(&st);
45  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
46  }
47  }
48 
49  void TransformationBase::set_matrix_mode(const STMatMode mode)
50  {
51  int ierr = STSetMatMode(st,mode);
52  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
53  }
54 
56  {
57  int 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  int ierr;
76  ierr = STSetType (st, const_cast<char *>(STSHIFT));
77  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
78 
79  ierr = STSetShift (st, additional_data.shift_parameter);
80  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
81  }
82 
83  /* ---------------- TransformationShiftInvert ------------------ */
84 
86  AdditionalData (const double shift_parameter)
87  :
88  shift_parameter (shift_parameter)
89  {}
90 
91  TransformationShiftInvert::TransformationShiftInvert (const MPI_Comm &mpi_communicator,
92  const AdditionalData &data)
93  :
94  TransformationBase(mpi_communicator),
95  additional_data (data)
96  {
97  int ierr;
98 #if DEAL_II_PETSC_VERSION_LT(3,1,0)
99  ierr = STSetType (st, const_cast<char *>(STSINV));
100 #else
101  ierr = STSetType (st, const_cast<char *>(STSINVERT));
102 #endif
103  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
104 
105  ierr = STSetShift (st, additional_data.shift_parameter);
106  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
107  }
108 
109  /* --------------- TransformationSpectrumFolding ----------------- */
110 
112  AdditionalData (const double shift_parameter)
113  :
114  shift_parameter (shift_parameter)
115  {}
116 
118  const AdditionalData &data)
119  :
120  TransformationBase(mpi_communicator),
121  additional_data (data)
122  {
123 #if DEAL_II_PETSC_VERSION_LT(3,5,0)
124  int ierr;
125  ierr = STSetType (st, const_cast<char *>(STFOLD));
126  (void)ierr;
127  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
128 
129  ierr = STSetShift (st, additional_data.shift_parameter);
130  (void)ierr;
131  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
132 #else
133  // PETSc/SLEPc version must be < 3.5.0.
134  (void)st;
135  Assert ((false),
136  ExcMessage ("Folding transformation has been removed in SLEPc 3.5.0 and newer."
137  "You cannot use this transformation anymore."));
138 #endif
139  }
140 
141  /* ------------------- TransformationCayley --------------------- */
142 
144  AdditionalData (const double shift_parameter,
145  const double antishift_parameter)
146  :
147  shift_parameter (shift_parameter),
148  antishift_parameter (antishift_parameter)
149  {
150  }
151 
152  TransformationCayley::TransformationCayley (const MPI_Comm &mpi_communicator,
153  const AdditionalData &data)
154  :
155  TransformationBase(mpi_communicator),
156  additional_data (data)
157  {
158  int ierr = STSetType (st, const_cast<char *>(STCAYLEY));
159  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
160 
161  ierr = STSetShift (st, additional_data.shift_parameter);
162  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
163 
164  ierr = STCayleySetAntishift (st, additional_data.antishift_parameter);
165  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
166  }
167 
168 }
169 
170 DEAL_II_NAMESPACE_CLOSE
171 
172 #endif // DEAL_II_WITH_SLEPC
void set_solver(const PETScWrappers::SolverBase &solver)
::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< SolverData > solver_data
Definition: petsc_solver.h:250
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
TransformationShiftInvert(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
#define Assert(cond, exc)
Definition: exceptions.h:294
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)
TransformationSpectrumFolding(const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())