dune-fem  2.4.1-rc
petsccommon.hh
Go to the documentation of this file.
1 // vim: set expandtab ts=2 sw=2 sts=2:
2 #ifndef DUNE_FEM_PETSCCOMMON_HH
3 #define DUNE_FEM_PETSCCOMMON_HH
4 
5 /*
6  * This file should be the first include in every dune-fem-petsc related header
7  * It contains some common code in order to deal with PETSc in a consistent, object oriented
8  * fashion
9  */
10 
11 #include <string>
12 #include <iostream>
13 
14 #include <dune/common/stdstreams.hh>
15 #include <dune/common/exceptions.hh>
16 
17 #if HAVE_PETSC
18 
19  /*
20  * This turns off PETSc logging of MPI calls. If it is on, PETSc redifines MPI functions as macros,
21  * which breaks some code. E.g. MPI_Allreduce is redefined
22  */
23 #define PETSC_HAVE_BROKEN_RECURSIVE_MACRO
24 
25 #include <petsc.h>
26 
27 namespace Dune
28 {
29 
30  namespace Petsc
31  {
32  // PETSc #defines VecType to be char* - this can cause problems with e.g. ALUGrid 1.50. So we #undef
33  // it here and use a typedef instead. The following trick is really nasty, but conforming to the standard and
34  // - most important of all - working :)
35  typedef VecType
36 #undef VecType
37  VecType;
38 
39 
40 #if HAVE_MPI
41 #define FEM_PETSC_COMM_DEFAULT PETSC_COMM_WORLD
42 #else
43 #define FEM_PETSC_COMM_DEFAULT PETSC_COMM_SELF
44 #endif
45 
46  /*
47  * exceptions
48  */
49  class Exception : public Dune::Exception {};
50 
51  /*
52  * types
53  */
54  typedef ::PetscErrorCode ErrorCode;
55 
56 
57 
58  /* The following 2 methods are needed to make the code work with the CHKERRQ
59  * macro (which we want to use because it can print useful diagnostic output in
60  * case of errors)
61  */
62  inline ErrorCode ErrorCheckHelper ( ErrorCode errorCode ) { CHKERRQ( errorCode ); return 0; }
63 
64  inline ErrorCode ErrorHandler ( MPI_Comm comm, int line, const char *function, const char *file,
65 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
66  const char *dir,
67 #endif
68  ErrorCode errorCode, PetscErrorType p, const char *message, void *context )
69  {
70  std::ostringstream msgout;
71  msgout << "PETSc Error in the PETSc function '" << function << "' at " << file << ":" << line << ":";
72 
73  const char *text;
74  PetscErrorMessage( errorCode, &text, PETSC_NULL );
75  if( text )
76  msgout << " '" << text << "'. Error message: '" << message << "'";
77  else
78  msgout << " Unable to retrieve error text";
79 
80  DUNE_THROW( Exception, msgout.str() );
81 
82  return errorCode;
83  }
84 
85  inline void ErrorCheck ( ErrorCode errorCode )
86  {
87  if( errorCode )
88  {
89  DUNE_THROW( Exception, "Generic PETSC exception" );
90  }
91  }
92 
93  /*
94  * This should be called right after the initialization of the MPI manager. It expects the same arguments
95  * as PetscInitialize
96  */
97  inline void initialize( const bool verbose, int &argc, char **&args, const char file[] = 0 , const char help[] = 0, bool ownHandler = true )
98  {
99  ::PetscInitialize( &argc, &args, file, help );
100 
101  if( ownHandler )
102  {
103  // set up our error handler
104  if( verbose )
105  {
106  dverb << "INFORMATION: Setting up an own error handler for PETSc errors. If you want the default PETSc handler,\n"
107  << "INFORMATION: set the last argument of Dune::Petsc::initialize(...) to 'false'.\n";
108  }
109  ::PetscPushErrorHandler( &ErrorHandler, 0 );
110  }
111  }
112 
113 #if 0
114  /*
115  * This should be called right after the initialization of the MPI manager. It expects the same arguments
116  * as PetscInitialize
117  */
118  inline void initialize( int *argc, char ***args, const char file[], const char help[], bool ownHandler = true )
119  {
120  ::PetscInitialize( argc, args, file, help );
121 
122  if( ownHandler )
123  {
124  // set up our error handler
125  if( Dune::Fem::MPIManager::rank() == 0 )
126  {
127  std::cerr << "INFORMATION: Setting up an own error handler for PETSc errors. If you want the default PETSc handler,\n"
128  << "INFORMATION: set the last argument of Dune::Petsc::initialize(...) to 'false'.\n";
129  }
130  ::PetscPushErrorHandler( &ErrorHandler, 0 );
131  }
132  }
133 #endif
134 
135  /*
136  * This should be called just before the termination of the program.
137  */
138  inline void finalize ()
139  {
140  // TODO: test here if we are using our own handler
141  ::PetscPopErrorHandler();
142 
143  ::PetscFinalize();
144  }
145 
146  /*
147  * ==================================================
148  * These are simple mappings to PETSc's C-routines. (Maybe some of them are not needed...).
149  *
150  * The PETSC_VERSION_... customizations are not very well tested yet
151  */
152  inline void KSPCreate ( KSP *inksp ) { ErrorCheck( ::KSPCreate( FEM_PETSC_COMM_DEFAULT, inksp ) ); }
153  inline void KSPDestroy ( KSP *ksp ) {
154 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
155  ErrorCheck( ::KSPDestroy( *ksp ) );
156 #else
157  ErrorCheck( ::KSPDestroy( ksp ) );
158 #endif // PETSC_PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
159  }
160  inline void KSPGetPC ( KSP ksp, PC *pc ) { ErrorCheck( ::KSPGetPC( ksp, pc ) ); }
161  inline void KSPSetFromOptions ( KSP ksp ) { ErrorCheck( ::KSPSetFromOptions( ksp ) ); }
162  inline void KSPSetType ( KSP ksp, const KSPType type ) { ErrorCheck( ::KSPSetType( ksp, type ) ); }
163  inline void KSPGMRESSetRestart ( KSP ksp, PetscInt restart ) { ErrorCheck( ::KSPGMRESSetRestart( ksp, restart ) ); }
164  inline void KSPView ( KSP ksp, PetscViewer viewer = PETSC_VIEWER_STDOUT_(FEM_PETSC_COMM_DEFAULT) )
165  {
166  ErrorCheck( ::KSPView( ksp, viewer ) );
167  }
168  inline void KSPMonitorSet (KSP ksp, PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),
169 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
170  void *mctx,PetscErrorCode (*monitordestroy)(void*)
171 #else
172  void *mctx,PetscErrorCode (*monitordestroy)(void**)
173 #endif
174  )
175  {
176  ErrorCheck( ::KSPMonitorSet( ksp, monitor, mctx, monitordestroy ) );
177  }
178  inline void KSPGetIterationNumber( KSP ksp, PetscInt* its )
179  { ErrorCheck( ::KSPGetIterationNumber( ksp, its ) ); }
180 
181 
182 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
183  inline void KSPSetOperators (KSP ksp, Mat Amat, Mat Pmat, MatStructure flag ) { ErrorCheck( ::KSPSetOperators( ksp, Amat, Pmat, flag ) ); }
184 #else
185  inline void KSPSetOperators (KSP ksp, Mat Amat, Mat Pmat ) { ErrorCheck( ::KSPSetOperators( ksp, Amat, Pmat ) ); }
186 #endif
187  inline void KSPSetTolerances ( KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits )
188  { ErrorCheck( ::KSPSetTolerances( ksp, rtol, abstol, dtol, maxits ) ); }
189  inline void KSPSolve ( KSP ksp, Vec b, Vec x ) { ErrorCheck( ::KSPSolve( ksp, b, x ) ); }
190  inline void KSPSetPC ( KSP ksp, PC pc ) { ErrorCheck( ::KSPSetPC( ksp, pc ) ); }
191 
192  // preconditioning
193  inline void PCCreate ( PC* pc) { ErrorCheck( ::PCCreate( FEM_PETSC_COMM_DEFAULT, pc ) ); }
194  inline void PCDestroy ( PC* pc) {
195 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
196  ErrorCheck( ::PCDestroy( *pc ) );
197 #else
198  ErrorCheck( ::PCDestroy( pc ) );
199 #endif
200  }
201  inline void PCSetType ( PC pc, const PCType type ) { ErrorCheck( ::PCSetType( pc, type ) ); }
202  inline void PCFactorSetLevels( PC pc, PetscInt level ) { ErrorCheck( ::PCFactorSetLevels( pc, level ) ); }
203  inline void PCFactorSetMatSolverPackage( PC pc, const MatSolverPackage type )
204  {
205  ErrorCheck( ::PCFactorSetMatSolverPackage( pc, type ) );
206  }
207  inline void PCHYPRESetType( PC pc, const char* type )
208  {
209 #if HAVE_PETSC_HYPRE
210  ErrorCheck( ::PCHYPRESetType( pc, type ) );
211 #else
212  DUNE_THROW( InvalidStateException, "HYPRE not found within PETSc, please re-configure!");
213 #endif
214  }
215 
216  // matrix routines
217  inline void MatAssemblyBegin ( Mat mat, MatAssemblyType type ) { ErrorCheck( ::MatAssemblyBegin( mat, type ) ); }
218  inline void MatAssemblyEnd ( Mat mat, MatAssemblyType type ) { ErrorCheck( ::MatAssemblyEnd( mat, type ) ); }
219  inline void MatCreate ( Mat *A ) { ErrorCheck( ::MatCreate( FEM_PETSC_COMM_DEFAULT, A) ); }
220  inline void MatCreateBlockMat ( Mat *A, PetscInt m, PetscInt n, PetscInt bs, PetscInt nz, PetscInt* nnz )
221  {
222  ErrorCheck( ::MatCreateBlockMat( FEM_PETSC_COMM_DEFAULT, n, m, bs, nz, nnz, A) );
223  }
224  inline void MatDestroy ( Mat *A ) {
225  #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
226  ErrorCheck( ::MatDestroy( *A ) );
227  #else
228  ErrorCheck( ::MatDestroy( A ) );
229  #endif // PETSC_PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
230  }
231  inline void MatSetUp( Mat mat )
232  {
233  ErrorCheck( ::MatSetUp(mat));
234  }
235  inline void MatSetUp( Mat mat, PetscInt bs, int nz )
236  {
237  if (bs == 1)
238  {
239  ErrorCheck( ::MatSeqAIJSetPreallocation(mat,nz,PETSC_NULL) );
240  ErrorCheck( ::MatMPIAIJSetPreallocation(mat,nz,PETSC_NULL,nz/2,PETSC_NULL) );
241  }
242  else
243  {
244  ErrorCheck( ::MatSeqBAIJSetPreallocation(mat,bs,nz,PETSC_NULL) );
245  ErrorCheck( ::MatMPIBAIJSetPreallocation(mat,bs,nz,PETSC_NULL,nz/2,PETSC_NULL) );
246  }
247  // the following seems not to work for block matrix
248  ErrorCheck( ::MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) );
249  // the following only works for block matrix
250  // but should be used with MAT_NEW_NONZERO_LOCATIONS...
251  // ErrorCheck( ::MatSetOption(mat, MAT_USE_HASH_TABLE,PETSC_FALSE) );
252  ErrorCheck( ::MatSetUp(mat));
253  }
254  inline void MatSetUp( Mat mat, PetscInt bs, const int *d_nnz, const int *o_nnz )
255  {
256  if (bs == 1)
257  {
258  ErrorCheck( ::MatSeqAIJSetPreallocation(mat,0,d_nnz ) );
259  ErrorCheck( ::MatMPIAIJSetPreallocation(mat,0,d_nnz,5,o_nnz) );
260  }
261  else
262  {
263  ErrorCheck( ::MatSeqBAIJSetPreallocation(mat,bs,0,d_nnz ) );
264  ErrorCheck( ::MatMPIBAIJSetPreallocation(mat,bs,0,d_nnz,5,PETSC_NULL) );
265  }
266  // see previous comments
267  ErrorCheck( ::MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) );
268  ErrorCheck( ::MatSetUp(mat));
269  }
270 
271  inline void MatGetOwnershipRange ( Mat mat, PetscInt *m, PetscInt* n ) { ErrorCheck( ::MatGetOwnershipRange( mat, m, n ) ); }
272  inline void MatGetSize ( Mat mat, PetscInt *m, PetscInt* n ) { ErrorCheck( ::MatGetSize( mat, m, n ) ); }
273  inline void MatMult ( Mat mat, Vec x, Vec y ) { ErrorCheck( ::MatMult( mat, x, y ) ); }
274  inline void MatSetBlockSize ( Mat A, PetscInt bs ) { ErrorCheck( ::MatSetBlockSize( A, bs ) ); }
275  inline void MatSetSizes ( Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N ) { ErrorCheck( ::MatSetSizes( A, m, n, M, N ) ); }
276  inline void MatSetFromOptions ( Mat B ) { ErrorCheck( ::MatSetFromOptions( B ) ); }
277  inline void MatSetType ( Mat mat, const MatType matype ) { ErrorCheck( ::MatSetType( mat, matype ) ); }
278  inline void MatSetValue ( Mat v, PetscInt i, PetscInt j, PetscScalar va, InsertMode mode ) { ErrorCheck( ::MatSetValue( v, i, j, va, mode ) ); }
279  inline void MatSetValues ( Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv )
280  { ErrorCheck( ::MatSetValues( mat, m, idxm, n, idxn, v, addv ) ); }
281  inline void MatGetValues ( Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[] )
282  { ErrorCheck( ::MatGetValues( mat, m, idxm, n, idxn, v ) ); }
283  inline void MatView ( Mat mat, PetscViewer viewer ) { ErrorCheck( ::MatView( mat, viewer ) ); }
284  inline void MatZeroEntries ( Mat mat ) { ErrorCheck( ::MatZeroEntries( mat ) ); }
285  inline void PetscBarrier ( PetscObject obj ) { ErrorCheck( ::PetscBarrier( obj ) ); }
286  inline void PetscFinalize () { ErrorCheck( ::PetscFinalize() ); }
287  inline void PetscInitialize( int *argc, char ***args, const char file[], const char help[] ) { ErrorCheck( ::PetscInitialize( argc, args, file, help ) ); }
288  inline void PetscViewerASCIIOpen ( MPI_Comm comm, const char name[], PetscViewer *lab ) { ErrorCheck( ::PetscViewerASCIIOpen( comm, name, lab ) ); }
289  inline void PetscViewerDestroy ( PetscViewer *viewer ) {
290  #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
291  ErrorCheck( ::PetscViewerDestroy( *viewer ) );
292  #else
293  ErrorCheck( ::PetscViewerDestroy( viewer ) );
294  #endif
295  }
296  inline void PetscViewerSetFormat ( PetscViewer viewer, PetscViewerFormat format ) { ErrorCheck( ::PetscViewerSetFormat( viewer, format ) ); }
297  inline void VecAssemblyBegin ( Vec vec ) { ErrorCheck( ::VecAssemblyBegin( vec ) ); }
298  inline void VecAssemblyEnd ( Vec vec ) { ErrorCheck( ::VecAssemblyEnd( vec ) ); }
299  inline void VecAXPY ( Vec y, PetscScalar alpha, Vec x) { ErrorCheck( ::VecAXPY( y, alpha, x ) ); }
300  inline void VecCopy ( Vec x, Vec y ) { ErrorCheck( ::VecCopy( x, y ) ); }
301  inline void VecCreate ( Vec *vec ) { ErrorCheck( ::VecCreate( FEM_PETSC_COMM_DEFAULT, vec ) ); }
302  inline void VecCreateGhost ( PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[], Vec *vv )
303  { ErrorCheck( ::VecCreateGhost( FEM_PETSC_COMM_DEFAULT, n, N, nghost, ghosts, vv ) ); }
304  inline void VecDestroy ( Vec *v ) {
305  #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
306  ErrorCheck( ::VecDestroy( *v ) );
307  #else
308  ErrorCheck( ::VecDestroy( v ) );
309  #endif // PETSC_PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
310  }
311  inline void VecDot ( Vec x, Vec y, PetscScalar *val ) { ErrorCheck( ::VecDot( x, y, val ) ); }
312  inline void VecDuplicate ( Vec v, Vec *newv ) { ErrorCheck( ::VecDuplicate( v, newv ) ); }
313  inline void VecGetLocalSize ( Vec x, PetscInt *size ) { ErrorCheck( ::VecGetLocalSize( x, size ) ); }
314  inline void VecGetOwnershipRange ( Vec x, PetscInt *low, PetscInt *high ) { ErrorCheck( ::VecGetOwnershipRange( x, low, high ) ); }
315  inline void VecGetSize ( Vec x, PetscInt *size ) { ErrorCheck( ::VecGetSize( x, size ) ); }
316  inline void VecGetValues ( Vec x, PetscInt ni, const PetscInt ix[], PetscScalar y[] ) { ErrorCheck( ::VecGetValues( x, ni, ix, y ) ); }
317  inline void VecGhostGetLocalForm ( Vec g, Vec *l ) { ErrorCheck( ::VecGhostGetLocalForm( g, l ) ); }
318  inline void VecGhostRestoreLocalForm ( Vec g, Vec *l ) { ErrorCheck( ::VecGhostRestoreLocalForm( g, l ) ); }
319  inline void VecGhostUpdateBegin ( Vec g, InsertMode insertmode, ScatterMode scattermode ) { ErrorCheck( ::VecGhostUpdateBegin( g, insertmode, scattermode ) ); }
320  inline void VecGhostUpdateEnd ( Vec g, InsertMode insertmode, ScatterMode scattermode ) { ErrorCheck( ::VecGhostUpdateEnd( g, insertmode, scattermode ) ); }
321  inline void VecNorm ( Vec x, NormType type, PetscReal *val ) { ErrorCheck( ::VecNorm( x, type, val ) ); }
322  inline void VecScale ( Vec x, PetscScalar alpha ) { ErrorCheck( ::VecScale( x, alpha ) ); }
323  inline void VecSet ( Vec x, PetscScalar alpha ) { ErrorCheck( ::VecSet( x, alpha ) ); }
324  inline void VecSetFromOptions ( Vec vec ) { ErrorCheck( ::VecSetFromOptions( vec ) ); }
325  inline void VecSetType ( Vec vec, const VecType method ) { ErrorCheck( ::VecSetType( vec, method ) ); }
326  inline void VecSetSizes ( Vec v, PetscInt n, PetscInt N ) { ErrorCheck( ::VecSetSizes( v, n, N ) ); }
327  inline void VecSetValue ( Vec v, int row, PetscScalar value, InsertMode mode ) { ErrorCheck( ::VecSetValue( v, row, value, mode ) ); }
328  inline void VecView ( Vec vec, PetscViewer viewer ) { ErrorCheck( ::VecView( vec, viewer ) ); }
329 
330  } // namespace Petsc
331 
332 } // namespace Dune
333 
334 #endif // #if HAVE_PETSC
335 
336 #endif // DUNE_FEM_PETSCCOMMON_HH
static int rank()
Definition: mpimanager.hh:116
int MPI_Comm
Definition: odeobjectfiles.cc:12
Definition: coordinate.hh:4