OOFEM 3.0
Loading...
Searching...
No Matches
slepcsolver.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35#include "slepcsolver.h"
36#include "petscsparsemtrx.h"
37#include "engngm.h"
38#include "floatarray.h"
39#include "verbose.h"
40#include "domain.h"
41#include "classfactory.h"
42
43#define TIME_REPORT
44
45#ifdef TIME_REPORT
46 #include "timer.h"
47#endif
48
49namespace oofem {
51
52
53SLEPcSolver :: SLEPcSolver(Domain *d, EngngModel *m) : SparseGeneralEigenValueSystemNM(d, m),
54 //A(nullptr),
55 //B(nullptr),
56 epsInit(false)
57{
58}
59
60
61SLEPcSolver :: ~SLEPcSolver()
62{
63 if ( epsInit ) {
64 EPSDestroy(& eps);
65 }
66}
67
69SLEPcSolver :: solve(SparseMtrx &a, SparseMtrx &b, FloatArray &_eigv, FloatMatrix &_r, double rtol, int nroot)
70{
71 PetscErrorCode ierr;
72 int size;
73 ST st;
74
75 // first check whether Lhs is defined
76
77 if ( a.giveNumberOfRows() != a.giveNumberOfColumns() ||
80 OOFEM_ERROR("matrices size mismatch");
81 }
82
83 PetscSparseMtrx* A = dynamic_cast< PetscSparseMtrx * >(&a);
84 PetscSparseMtrx* B = dynamic_cast< PetscSparseMtrx * >(&b);
85
86 if ( !A || !B ) {
87 OOFEM_ERROR("PetscSparseMtrx Expected");
88 }
89
90 size = engngModel->giveParallelContext( A->giveDomainIndex() )->giveNumberOfNaturalEqs(); // A->giveLeqs();
91
92 _r.resize(size, nroot);
93 _eigv.resize(nroot);
94
95
96 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97 * Create the eigensolver and set various options
98 * - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99 int nconv, nite;
100 EPSConvergedReason reason;
101
102#ifdef TIME_REPORT
103 Timer timer;
104 timer.startTimer();
105#endif
106
107 if ( !epsInit ) {
108 /*
109 * Create eigensolver context
110 */
111#ifdef __MPI_PARALLEL_MODE
112 MPI_Comm comm = engngModel->giveParallelComm();
113#else
114 MPI_Comm comm = PETSC_COMM_SELF;
115#endif
116 ierr = EPSCreate(comm, & eps);
117 CHKERRQ(ierr);
118 epsInit = true;
119 }
120
121 /*
122 * Set operators. In this case, it is a generalized eigenvalue problem
123 */
124
125 ierr = EPSSetOperators( eps, * A->giveMtrx(), * B->giveMtrx() );
126 CHKERRQ(ierr);
127 ierr = EPSSetProblemType(eps, EPS_GHEP);
128 CHKERRQ(ierr);
129 ierr = EPSGetST(eps, & st);
130 CHKERRQ(ierr);
131 ierr = STSetType(st, STSINVERT);
132 CHKERRQ(ierr);
133 ierr = STSetMatStructure(st, SAME_NONZERO_PATTERN);
134 CHKERRQ(ierr);
135 ierr = EPSSetTolerances(eps, ( PetscReal ) rtol, PETSC_DECIDE);
136 CHKERRQ(ierr);
137 ierr = EPSSetDimensions(eps, ( PetscInt ) nroot, PETSC_DECIDE, PETSC_DECIDE);
138 CHKERRQ(ierr);
139 ierr = EPSSetWhichEigenpairs(eps, EPS_SMALLEST_MAGNITUDE);
140 CHKERRQ(ierr);
141
142 /*
143 * Set solver parameters at runtime
144 */
145
146 ierr = EPSSetFromOptions(eps);
147 CHKERRQ(ierr);
148
149 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150 * Solve the eigensystem
151 * - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152
153 ierr = EPSSolve(eps);
154 CHKERRQ(ierr);
155
156 ierr = EPSGetConvergedReason(eps, & reason);
157 CHKERRQ(ierr);
158 ierr = EPSGetIterationNumber(eps, & nite);
159 CHKERRQ(ierr);
160 OOFEM_LOG_INFO("SLEPcSolver::solve EPSConvergedReason: %d, number of iterations: %d\n", reason, nite);
161
162 ierr = EPSGetConverged(eps, & nconv);
163 CHKERRQ(ierr);
164
165 if ( nconv > 0 ) {
166 OOFEM_LOG_INFO("SLEPcSolver :: solveYourselfAt: Convergence reached for RTOL=%20.15f", rtol);
167 PetscScalar kr;
168 Vec Vr;
169
170 ierr = MatGetVecs(* B->giveMtrx(), PETSC_NULL, & Vr);
171 CHKERRQ(ierr);
172
173 FloatArray Vr_loc;
174
175 for ( int i = 0; i < nconv && i < nroot; i++ ) {
176 ierr = EPSGetEigenpair(eps, i, & kr, PETSC_NULL, Vr, PETSC_NULL);
177 CHKERRQ(ierr);
178
179 //Store the eigenvalue
180 _eigv.at(i + 1) = kr;
181
182 //Store the eigenvector
183 A->scatterG2L(Vr, Vr_loc);
184 for ( int j = 0; j < size; j++ ) {
185 _r.at(j + 1, i + 1) = Vr_loc.at(j + 1);
186 }
187 }
188
189 ierr = VecDestroy(& Vr);
190 CHKERRQ(ierr);
191 } else {
192 OOFEM_ERROR("No converged eigenpairs");
193 }
194
195#ifdef TIME_REPORT
196 timer.stopTimer();
197 OOFEM_LOG_INFO( "SLEPcSolver info: user time consumed by solution: %.2fs\n", timer.getUtime() );
198#endif
199
200 return CR_CONVERGED;
201}
202} // end namespace oofem
#define REGISTER_GeneralizedEigenValueSolver(class, type)
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
double at(std::size_t i, std::size_t j) const
EngngModel * engngModel
Pointer to engineering model.
Definition nummet.h:86
int scatterG2L(Vec src, FloatArray &dest) const
Scatters global vector to natural one.
bool epsInit
Flag if context initialized.
Definition slepcsolver.h:61
EPS eps
Eigenvalue solver context.
Definition slepcsolver.h:59
SparseGeneralEigenValueSystemNM(Domain *d, EngngModel *m)
Constructor.
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition sparsemtrx.h:114
int giveNumberOfColumns() const
Returns number of columns of receiver.
Definition sparsemtrx.h:116
double getUtime()
Returns total user time elapsed in seconds.
Definition timer.C:105
void startTimer()
Definition timer.C:69
void stopTimer()
Definition timer.C:77
#define OOFEM_ERROR(...)
Definition error.h:79
#define OOFEM_LOG_INFO(...)
Definition logger.h:143

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011