OOFEM 3.0
Loading...
Searching...
No Matches
petscsolver.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 "petscsolver.h"
36
37#include "petscsparsemtrx.h"
38#include "convergedreason.h"
39#include "engngm.h"
40#include "floatarray.h"
41#include "verbose.h"
42#include "timer.h"
43#include "error.h"
44#include "classfactory.h"
45
46#include <petscksp.h>
47
48namespace oofem {
50
51PetscSolver :: PetscSolver(Domain *d, EngngModel *m) : SparseLinearSystemNM(d, m) { }
52
53PetscSolver :: ~PetscSolver() { }
54
55ConvergedReason PetscSolver :: solve(SparseMtrx &A, FloatArray &b, FloatArray &x)
56{
57 int neqs = b.giveSize();
58 if ( x.giveSize() != neqs )
59 x.resize(neqs);
60
61 PetscSparseMtrx *Lhs = dynamic_cast< PetscSparseMtrx * >(&A);
62 if ( !Lhs ) {
63 OOFEM_ERROR("PetscSparseMtrx Expected");
64 }
65
66 Vec globRhsVec;
67 Vec globSolVec;
68
69 /*
70 * scatter and gather rhs to global representation (automatically detects sequential/parallel modes)
71 */
72 Lhs->createVecGlobal(& globRhsVec);
73 Lhs->scatterL2G(b, globRhsVec);
74
75 Lhs->createVecGlobal(& globSolVec);
76 Lhs->scatterL2G(x, globSolVec);
77
78 //VecView(globRhsVec,PETSC_VIEWER_STDOUT_WORLD);
79 ConvergedReason s = this->petsc_solve(*Lhs, globRhsVec, globSolVec);
80 //VecView(globSolVec,PETSC_VIEWER_STDOUT_WORLD);
81
82 Lhs->scatterG2L(globSolVec, x);
83
84 VecDestroy(& globSolVec);
85 VecDestroy(& globRhsVec);
86
87 return s;
88}
89
91PetscSolver :: petsc_solve(PetscSparseMtrx &Lhs, Vec b, Vec x)
92{
93 int nite;
94 PetscErrorCode err;
95 KSPConvergedReason reason;
96
97 Timer timer;
98 timer.startTimer();
99
100 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101 * Create the linear solver and set various options
102 * - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103#ifdef _OPENMP
104 #pragma omp critical
105#endif
106 if ( !Lhs.kspInit ) {
107#ifdef __MPI_PARALLEL_MODE
108 MPI_Comm comm = engngModel->giveParallelComm();
109#else
110 MPI_Comm comm = PETSC_COMM_SELF;
111#endif
112 KSPCreate(comm, & Lhs.ksp);
113 Lhs.kspInit = true;
114 }
115
116 /*
117 * Set operators. Here the matrix that defines the linear system
118 * also serves as the preconditioning matrix.
119 */
120 if ( Lhs.newValues ) { // Optimization for successive solves
122 if ( this->engngModel->requiresUnknownsDictionaryUpdate() ) {
123#if PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR >= 5
124 // The syntax for this function changed in PETSc version 3.5.1 /ES
125 KSPSetOperators(Lhs.ksp, * Lhs.giveMtrx(), * Lhs.giveMtrx());
126#else
127 KSPSetOperators(Lhs.ksp, * Lhs.giveMtrx(), * Lhs.giveMtrx(), DIFFERENT_NONZERO_PATTERN);
128#endif
129 } else {
130 //KSPSetOperators(Lhs->ksp, * Lhs->giveMtrx(), * Lhs->giveMtrx(), SAME_PRECONDITIONER);
131#if PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR >= 5
132 // The syntax for this function changed in PETSc version 3.5.1 /ES
133 KSPSetOperators(Lhs.ksp, * Lhs.giveMtrx(), * Lhs.giveMtrx());
134#else
135 KSPSetOperators(Lhs.ksp, * Lhs.giveMtrx(), * Lhs.giveMtrx(), SAME_NONZERO_PATTERN);
136#endif
137 }
138 Lhs.newValues = false;
139
140 /*
141 * Set linear solver defaults for this problem (optional).
142 * - The following two statements are optional; all of these
143 * parameters could alternatively be specified at runtime via
144 * KSPSetFromOptions(). All of these defaults can be
145 * overridden at runtime, as indicated below.
146 */
147 KSPSetTolerances(Lhs.ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
148
149 /*
150 * Set runtime options, e.g.,
151 * -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
152 * These options will override those specified above as long as
153 * KSPSetFromOptions() is called _after_ any other customization
154 * routines.
155 */
156 KSPSetFromOptions(Lhs.ksp);
157 }
158
159 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160 * Solve the linear system
161 * - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162 //MatView(*Lhs->giveMtrx(),PETSC_VIEWER_STDOUT_SELF);
163 err = KSPSolve(Lhs.ksp, b, x);
164 if ( err != 0 ) {
165 OOFEM_ERROR("Error when solving: %d\n", err);
166 }
167 KSPGetConvergedReason(Lhs.ksp, & reason);
168 KSPGetIterationNumber(Lhs.ksp, & nite);
169
170 if ( reason >= 0 ) {
171 //OOFEM_LOG_INFO("PetscSolver: Converged. KSPConvergedReason: %d, number of iterations: %d\n", reason, nite);
172 } else {
173 OOFEM_WARNING("Diverged! KSPConvergedReason: %d, number of iterations: %d\n", reason, nite);
174 }
175
176 timer.stopTimer();
177 if ( this->engngModel->giveProblemScale() == macroScale ) {
178 OOFEM_LOG_INFO( "PetscSolver: User time consumed by solution: %.2fs, KSPConvergedReason: %d, number of iterations: %d\n", timer.getUtime(), reason, nite );
179 }
180
181 if ( reason < 0 ) {
182 return CR_DIVERGED_ITS;
183 } else {
184 return CR_CONVERGED;
185 }
186}
187
188
189#if 0
191ConvergedReason PetscSolver :: solve(SparseMtrx &A, FloatMatrix &B, FloatMatrix &X)
192{
193 PetscSparseMtrx *Lhs = dynamic_cast< PetscSparseMtrx * >(&A);
194 if ( !Lhs ) {
195 OOFEM_ERROR("PetscSparseMtrx Expected");
196 }
197
198 Vec globRhsVec;
199 Vec globSolVec;
200
201 bool newLhs = true;
202 int rows = B.giveNumberOfRows();
203 int cols = B.giveNumberOfColumns();
205 X.resize(rows, cols);
206 double *Xptr = X.givePointer();
207
208 for ( int i = 0; i < cols; ++i ) {
209 VecCreateSeqWithArray(PETSC_COMM_SELF, rows, B.givePointer() + rows * i, & globRhsVec);
210 VecDuplicate(globRhsVec, & globSolVec);
211 s = this->petsc_solve(*Lhs, globRhsVec, globSolVec, newLhs);
212 if ( !( s == CR_CONVERGED ) ) {
213 OOFEM_WARNING("No success at solving column %d", i + 1);
214 return s;
215 }
216 newLhs = false;
217 double *ptr;
218 VecGetArray(globSolVec, & ptr);
219 for ( int j = 0; j < rows; ++j ) {
220 Xptr [ j + rows * i ] = ptr [ j ];
221 }
222 VecRestoreArray(globSolVec, & ptr);
223 }
224 VecDestroy(& globSolVec);
225 VecDestroy(& globRhsVec);
226 return s;
227}
228#endif
229} // end namespace oofem
#define REGISTER_SparseLinSolver(class, type)
void resize(Index s)
Definition floatarray.C:94
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
const double * givePointer() const
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
int giveNumberOfColumns() const
Returns number of columns of receiver.
int giveNumberOfRows() const
Returns number of rows of receiver.
EngngModel * engngModel
Pointer to engineering model.
Definition nummet.h:86
ConvergedReason petsc_solve(PetscSparseMtrx &A, Vec b, Vec x)
Definition petscsolver.C:91
int scatterL2G(const FloatArray &src, Vec dest) const
bool kspInit
Flag if context initialized.
void createVecGlobal(Vec *answer) const
Creates a global vector that fits the instance of this matrix.
KSP ksp
Linear solver context.
bool newValues
Flag if matrix has changed since last solve.
int scatterG2L(Vec src, FloatArray &dest) const
Scatters global vector to natural one.
SparseLinearSystemNM(Domain *d, EngngModel *m)
Constructor.
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_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
@ CR_DIVERGED_ITS
@ macroScale
Definition problemmode.h:46

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