File: PETScLUSolver.cpp

package info (click to toggle)
dolfin 2018.1.0.post1-16
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 28,764 kB
  • sloc: xml: 104,040; cpp: 98,856; python: 22,511; makefile: 204; sh: 182
file content (384 lines) | stat: -rw-r--r-- 12,787 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
// Copyright (C) 2005-2017 Anders Logg and Garth N. Wells
//
// This file is part of DOLFIN.
//
// DOLFIN is free software: you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// DOLFIN is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.

#ifdef HAS_PETSC

#include <petscksp.h>
#include <petscpc.h>
#include <dolfin/common/constants.h>
#include <dolfin/common/MPI.h>
#include <dolfin/common/NoDeleter.h>
#include <dolfin/common/Timer.h>
#include <dolfin/log/log.h>
#include <dolfin/parameter/GlobalParameters.h>
#include "LUSolver.h"
#include "PETScMatrix.h"
#include "PETScObject.h"
#include "PETScVector.h"
#include "PETScLUSolver.h"

using namespace dolfin;

// Functions in anonymous namespace (local scope)
namespace
{
  const MatSolverType get_solver_package_type(KSP ksp)
  {
    PetscErrorCode ierr;
#if PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR <= 8 && PETSC_VERSION_RELEASE == 1
    const MatSolverType solver_type;
#else
    MatSolverType solver_type = 0;
#endif
    PC pc;
    ierr = KSPGetPC(ksp, &pc);
    if (ierr != 0) dolfin::PETScObject::petsc_error(ierr, __FILE__, "KSPGetPC");

    ierr = PCFactorGetMatSolverType(pc, &solver_type);
    if (ierr != 0) dolfin::PETScObject::petsc_error(ierr, __FILE__, "PCFactorGetMatSolverType");

    return solver_type;
  }
  //---------------------------------------------------------------------------
  std::map<const MatSolverType, bool> methods_cholesky
  = { {MATSOLVERUMFPACK,      false},
      {MATSOLVERMUMPS,        true},
      {MATSOLVERPASTIX,       true},
      {MATSOLVERSUPERLU,      false},
      {MATSOLVERSUPERLU_DIST, false},
      {MATSOLVERPETSC,        true} };

  //---------------------------------------------------------------------------
  bool solver_has_cholesky(const MatSolverType package)
  {
    auto it = methods_cholesky.find(package);
    dolfin_assert(it != methods_cholesky.end());
    return it->second;
  }
  //---------------------------------------------------------------------------
  const std::map<std::string, std::string> methods_descr
  = { {"default", "default LU solver"},
#if PETSC_HAVE_UMFPACK || PETSC_HAVE_SUITESPARSE
      {"umfpack", "UMFPACK (Unsymmetric MultiFrontal sparse LU factorization)"},
#endif
#if PETSC_HAVE_MUMPS
      {"mumps", "MUMPS (MUltifrontal Massively Parallel Sparse direct Solver)"},
#endif
#if PETSC_HAVE_PASTIX
      {"pastix", "PaStiX (Parallel Sparse matriX package)"},
#endif
#if PETSC_HAVE_SUPERLU
      {"superlu", "SuperLU"},
#endif
#if PETSC_HAVE_SUPERLU_DIST
      {"superlu_dist", "Parallel SuperLU"},
#endif
      {"petsc", "PETSc built in LU solver"} };

}
//-----------------------------------------------------------------------------

// List of available LU solvers
std::map<std::string, const MatSolverType> PETScLUSolver::lumethods
= { {"default", ""},
#if PETSC_HAVE_UMFPACK || PETSC_HAVE_SUITESPARSE
    {"umfpack",      MATSOLVERUMFPACK},
#endif
#if PETSC_HAVE_MUMPS
    {"mumps",        MATSOLVERMUMPS},
#endif
#if PETSC_HAVE_PASTIX
    {"pastix",       MATSOLVERPASTIX},
#endif
#if PETSC_HAVE_SUPERLU
      {"superlu",      MATSOLVERSUPERLU},
#endif
#if PETSC_HAVE_SUPERLU_DIST
    {"superlu_dist", MATSOLVERSUPERLU_DIST},
#endif
    {"petsc",        MATSOLVERPETSC}};

//-----------------------------------------------------------------------------
std::map<std::string, std::string> PETScLUSolver::methods()
{
  return methods_descr;
}
//-----------------------------------------------------------------------------
Parameters PETScLUSolver::default_parameters()
{
  Parameters p(LUSolver::default_parameters());
  p.rename("petsc_lu_solver");

  return p;
}
//-----------------------------------------------------------------------------
PETScLUSolver::PETScLUSolver(MPI_Comm comm, std::string method)
  :  PETScLUSolver(comm, nullptr, method)
{
  // Do nothing
}
//-----------------------------------------------------------------------------
PETScLUSolver::PETScLUSolver(std::string method)
  : PETScLUSolver(MPI_COMM_WORLD, nullptr, method)
{
  // Do nothing
}
//-----------------------------------------------------------------------------
PETScLUSolver::PETScLUSolver(MPI_Comm comm,
                             std::shared_ptr<const PETScMatrix> A,
                             std::string method) : _solver(comm)
{
  PetscErrorCode ierr;

  // Check dimensions, and check for symmetry
  PetscBool is_symmetric = PETSC_FALSE;
  if (A)
  {
    if (A->size(0) != A->size(1))
    {
      dolfin_error("PETScLUSolver.cpp",
                   "create PETSc LU solver",
                   "Cannot LU factorize non-square PETSc matrix");
    }

    dolfin_assert(A->mat());
    PetscBool symm_is_set = PETSC_FALSE;
    ierr = MatIsSymmetricKnown(A->mat(), &symm_is_set, &is_symmetric);
    if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "MatIsSymmetricKnown");
  }

  // Set parameter values
  parameters = default_parameters();

  // Select solver package
  const MatSolverType solver_package = select_solver(comm, method);

  // Get KSP pointer
  KSP ksp = _solver.ksp();

  // Make solver preconditioner only
  ierr = KSPSetType(ksp, KSPPREONLY);
  if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "KSPSetType");

  // Get PC
  PC pc;
  ierr = KSPGetPC(ksp, &pc);
  if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "KSPGetPC");

  // Set PC type to LU or PCCHOLESKY (depending on matrix symmetry)
  if (is_symmetric == PETSC_TRUE and solver_has_cholesky(solver_package))
  {
    ierr = PCSetType(pc, PCCHOLESKY);
    if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "PCSetType");
  }
  else
  {
    ierr = PCSetType(pc, PCLU);
    if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "PCSetType");
  }

  // Set LU solver package
  ierr = PCFactorSetMatSolverType(pc, solver_package);
  if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "PCFactorSetMatSolverType");

  // Set operator
  if (A)
  {
    ierr = KSPSetOperators(ksp, A->mat(), A->mat());
    if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "KSPSetOperators");
  }
}
//-----------------------------------------------------------------------------
PETScLUSolver::PETScLUSolver(std::shared_ptr<const PETScMatrix> A,
                             std::string method)
  : PETScLUSolver(MPI_COMM_WORLD, A, method)
{
  // Do nothing
}
//-----------------------------------------------------------------------------
PETScLUSolver::~PETScLUSolver()
{
  // Do nothing
}
//-----------------------------------------------------------------------------
void
PETScLUSolver::set_operator(std::shared_ptr<const GenericLinearOperator> A)
{
  _solver.set_operator(A);
}
//-----------------------------------------------------------------------------
void PETScLUSolver::set_operator(const PETScMatrix& A)
{
  _solver.set_operator(A);
}
//-----------------------------------------------------------------------------
std::size_t PETScLUSolver::solve(GenericVector& x, const GenericVector& b)
{
  return solve(x, b, false);
}
//-----------------------------------------------------------------------------
std::size_t PETScLUSolver::solve(GenericVector& x, const GenericVector& b,
                                 bool transpose)
{
  // FIXME: This should really go in PETScKrylovSolver

  const bool report = parameters["report"].is_set() ? parameters["report"] : false;
  if (report && dolfin::MPI::rank(mpi_comm()) == 0)
  {
    // Get PETSc operators
    Mat _A, _P;
    KSPGetOperators(_solver.ksp(), &_A, &_P);
    dolfin_assert(_A);
    PETScBaseMatrix A(_A);

    const MatSolverType solver_type = get_solver_package_type(_solver.ksp());
    log(PROGRESS,"Solving linear system of size %ld x %ld (PETSc LU solver, %s).",
        A.size(0), A.size(1), solver_type);
  }

  return _solver.solve(x, b);
}
//-----------------------------------------------------------------------------
std::size_t PETScLUSolver::solve(const GenericLinearOperator& A,
                                 GenericVector& x,
                                 const GenericVector& b)
{
  return solve(as_type<const PETScMatrix>(require_matrix(A)),
               as_type<PETScVector>(x),
               as_type<const PETScVector>(b));
}
//-----------------------------------------------------------------------------
std::size_t PETScLUSolver::solve(const PETScMatrix& A, PETScVector& x,
                                 const PETScVector& b)
{
  _solver.set_operators(A, A);
  return solve(x, b);
}
//-----------------------------------------------------------------------------
void PETScLUSolver::set_options_prefix(std::string options_prefix)
{
  _solver.set_options_prefix(options_prefix);
}
//-----------------------------------------------------------------------------
std::string PETScLUSolver::get_options_prefix() const
{
  return _solver.get_options_prefix();
}
//-----------------------------------------------------------------------------
void PETScLUSolver::set_from_options() const
{
  _solver.set_from_options();
}
//-----------------------------------------------------------------------------
MPI_Comm PETScLUSolver::mpi_comm() const
{
  return _solver.mpi_comm();
}
//-----------------------------------------------------------------------------
std::string PETScLUSolver::str(bool verbose) const
{
  std::stringstream s;

  if (verbose)
  {
    warning("Verbose output for PETScLUSolver not implemented, calling PETSc KSPView directly.");
    PetscErrorCode ierr = KSPView(_solver.ksp(), PETSC_VIEWER_STDOUT_WORLD);
    if (ierr != 0) PETScObject::petsc_error(ierr, __FILE__, "KSPView");
  }
  else
    s << "<PETScLUSolver>";

  return s.str();
}
//-----------------------------------------------------------------------------
KSP PETScLUSolver::ksp() const
{
  return _solver.ksp();
}
//-----------------------------------------------------------------------------
const MatSolverType PETScLUSolver::select_solver(MPI_Comm comm,
                                                 std::string method)
{
  // Check package string
  if (lumethods.count(method) == 0)
  {
    dolfin_error("PETScLUSolver.cpp",
                 "solve linear system using PETSc LU solver",
                 "Unknown LU method \"%s\"", method.c_str());
  }

  // Choose appropriate 'default' solver
  if (method == "default")
  {
    #if defined(PETSC_USE_64BIT_INDICES)
    if (dolfin::MPI::size(comm) == 1)
    {
      #if PETSC_HAVE_UMFPACK || PETSC_HAVE_SUITESPARSE
      method = "umfpack";
      // superlu_dist does not work efficiently since mc64 cannot be redistributed.  So use PETSc in preference to superlu_dist
      #else
      method = "petsc";
      warning("Using PETSc native LU solver. Consider configuring PETSc with an efficient LU solver (e.g. Umfpack).");
      #endif
    }
    else
    {
      // superlu_dist does not work efficiently since mc64 cannot be redistributed.  So use PETSc in preference to superlu_dist
      method = "petsc";
      warning("Using PETSc native LU solver. Consider specifying a more efficient LU solver (e.g.umfpack) if available.");
    }
    #else
    if (dolfin::MPI::size(comm) == 1)
    {
      #if PETSC_HAVE_UMFPACK || PETSC_HAVE_SUITESPARSE
      method = "umfpack";
      #elif PETSC_HAVE_MUMPS
      method = "mumps";
      #elif PETSC_HAVE_PASTIX
      method = "pastix";
      #elif PETSC_HAVE_SUPERLU
      method = "superlu";
      // superlu_dist does not work efficiently since mc64 cannot be redistributed.  So use PETSc in preference to superlu_dist
      #else
      method = "petsc";
      warning("Using PETSc native LU solver. Consider configuring PETSc with an efficient LU solver (e.g. UMFPACK, MUMPS).");
      #endif
    }
    else
    {
      #if PETSC_HAVE_MUMPS
      method = "mumps";
      #elif PETSC_HAVE_PASTIX
      method = "pastix";
      #elif PETSC_HAVE_SUPERLU_DIST
      method = "superlu_dist";
      #else
      dolfin_error("PETScLUSolver.cpp",
                   "solve linear system using PETSc LU solver",
                   "No suitable solver for parallel LU found. Consider configuring PETSc with MUMPS or SuperLU_dist");
#endif
    }
    #endif
  }

  auto it = lumethods.find(method);
  dolfin_assert(it !=  lumethods.end());
  return it->second;
}
//-----------------------------------------------------------------------------

#endif