DUNE-FEM (unstable)

petsccommon.hh
1#ifndef DUNE_FEM_PETSCCOMMON_HH
2#define DUNE_FEM_PETSCCOMMON_HH
3
4/*
5 * This file should be the first include in every dune-fem-petsc related header
6 * It contains some common code in order to deal with PETSc in a consistent, object oriented
7 * fashion
8 */
9
10#include <string>
11#include <iostream>
12#include <iomanip>
13
17
18#if HAVE_PETSC
19
20 /*
21 * This turns off PETSc logging of MPI calls. If it is on, PETSc redifines MPI functions as macros,
22 * which breaks some code. E.g. MPI_Allreduce is redefined
23 */
24#define PETSC_HAVE_BROKEN_RECURSIVE_MACRO
25
26#include <petsc.h>
27
28// newer versions of PETSc use PETSC_NULLPTR instread of PETSC_NULL
29#ifndef PETSC_NULLPTR
30#define PETSC_NULLPTR PETSC_NULL
31#endif
32
33namespace Dune
34{
35
36 namespace Petsc
37 {
38 // PETSc #defines VecType to be char* - this can cause problems with e.g. ALUGrid 1.50. So we #undef
39 // it here and use a typedef instead. The following trick is really nasty, but conforming to the standard and
40 // - most important of all - working :)
41 typedef VecType
42#undef VecType
43 VecType;
44
45 /*
46 * exceptions
47 */
48 class Exception : public Dune::Exception {};
49
50 /*
51 * types
52 */
53 typedef ::PetscErrorCode ErrorCode;
54
55
56
57 /* The following 2 methods are needed to make the code work with the CHKERRQ
58 * macro (which we want to use because it can print useful diagnostic output in
59 * case of errors)
60 */
61 inline static ErrorCode ErrorCheckHelper ( ErrorCode errorCode ) { CHKERRQ( errorCode ); return 0; }
62
63 inline ErrorCode ErrorHandler ( MPI_Comm comm, int line, const char *function, const char *file,
64#if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
65 const char *dir,
66#endif
67 ErrorCode errorCode, PetscErrorType p, const char *message, void *context )
68 {
69 std::ostringstream msgout;
70 msgout << "PETSc Error in the PETSc function '" << function << "' at " << file << ":" << line << ":";
71
72 const char *text;
73 PetscErrorMessage( errorCode, &text, PETSC_NULLPTR );
74 if( text )
75 msgout << " '" << text << "'. Error message: '" << message << "'";
76 else
77 msgout << " Unable to retrieve error text";
78
79 DUNE_THROW( Exception, msgout.str() );
80
81 return errorCode;
82 }
83
84 inline static void ErrorCheck ( ErrorCode errorCode )
85 {
86 if( errorCode )
87 {
88 DUNE_THROW( Exception, "Generic PETSC exception" );
89 }
90 }
91
92 /*
93 * This should be called right after the initialization of the MPI manager. It expects the same arguments
94 * as PetscInitialize
95 */
96 inline bool initialize( const bool verbose, int &argc, char **&args, const char file[] = 0 , const char help[] = 0, bool ownHandler = true )
97 {
98 bool wasInitializedHere = false ;
99 PetscBool petscInitialized = PETSC_FALSE;
100
101 // check whether PETSc had been initialized elsewhere
102 ::PetscInitialized( &petscInitialized );
103
104 if( ! petscInitialized )
105 {
106 ::PetscInitialize( &argc, &args, file, help );
107 wasInitializedHere = true;
108 }
109
110 if( ownHandler )
111 {
112 // set up our error handler
113 if( verbose )
114 {
115 dverb << "INFORMATION: Setting up an own error handler for PETSc errors. If you want the default PETSc handler,\n"
116 << "INFORMATION: set the last argument of Dune::Petsc::initialize(...) to 'false'.\n";
117 }
118 ::PetscPushErrorHandler( &ErrorHandler, 0 );
119 }
120 return wasInitializedHere;
121 }
122
123 /*
124 * This should be called just before the termination of the program.
125 */
126 inline void finalize ()
127 {
128 // TODO: test here if we are using our own handler
129 ::PetscPopErrorHandler();
130 PetscBool finalized = PETSC_FALSE;
131 ErrorCheck( ::PetscFinalized( &finalized ) );
132
133 if( ! finalized )
134 {
135 ::PetscFinalize();
136 }
137 }
138
139 template <class Comm>
140 inline auto castToPetscComm( const Comm& comm )
141 {
142 // this is needed because Dune::No_Comm or
143 // Dune::Communication< No_Comm > does not cast into MPI_Comm
144 if constexpr ( std::is_same< Dune::No_Comm, Comm > :: value ||
145 std::is_same< Dune::Communication< No_Comm >, Comm >::value )
146 {
147 return PETSC_COMM_SELF;
148 }
149 else
150 {
151 return comm;
152 }
153 }
154
155 /*
156 * ==================================================
157 * These are simple mappings to PETSc's C-routines. (Maybe some of them are not needed...).
158 *
159 * The PETSC_VERSION_... customizations are not very well tested yet
160 */
161 template <class Comm>
162 inline static void KSPCreate ( const Comm& comm, KSP *inksp )
163 {
164 ErrorCheck( ::KSPCreate( castToPetscComm( comm ), inksp ) );
165 }
166
167 inline static void KSPDestroy ( KSP *ksp ) {
168#if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
169 ErrorCheck( ::KSPDestroy( *ksp ) );
170#else
171 ErrorCheck( ::KSPDestroy( ksp ) );
172#endif // PETSC_PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
173 }
174 inline static void KSPGetPC ( KSP ksp, PC *pc ) { ErrorCheck( ::KSPGetPC( ksp, pc ) ); }
175 inline static void KSPSetFromOptions ( KSP ksp ) { ErrorCheck( ::KSPSetFromOptions( ksp ) ); }
176 inline static void KSPSetUp ( KSP ksp ) { ErrorCheck( ::KSPSetUp( ksp ) ); }
177 inline static void KSPSetType ( KSP ksp, const KSPType type ) { ErrorCheck( ::KSPSetType( ksp, type ) ); }
178 inline static void KSPGMRESSetRestart ( KSP ksp, PetscInt restart ) { ErrorCheck( ::KSPGMRESSetRestart( ksp, restart ) ); }
179
180 template <class Comm>
181 inline static void KSPView ( const Comm& comm,
182 KSP ksp,
183 PetscViewer viewer )
184 {
185 ErrorCheck( ::KSPView( ksp, viewer ) );
186 }
187
188 template <class Comm>
189 inline static void KSPView ( const Comm& comm,
190 KSP ksp )
191 {
192 KSPView( comm, ksp, PETSC_VIEWER_STDOUT_( castToPetscComm(comm) ) );
193 }
194
195 inline static void KSPMonitorSet (KSP ksp, PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),
196#if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
197 void *mctx,PetscErrorCode (*monitordestroy)(void*)
198#else
199 void *mctx,PetscErrorCode (*monitordestroy)(void**)
200#endif
201 )
202 {
203 ErrorCheck( ::KSPMonitorSet( ksp, monitor, mctx, monitordestroy ) );
204 }
205 inline static void KSPGetIterationNumber( KSP ksp, PetscInt* its )
206 { ErrorCheck( ::KSPGetIterationNumber( ksp, its ) ); }
207 inline static void KSPGetConvergedReason(KSP ksp, KSPConvergedReason *reason)
208 { ErrorCheck( ::KSPGetConvergedReason( ksp, reason ) ); }
209
210
211#if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
212 inline static void KSPSetOperators (KSP ksp, Mat Amat, Mat Pmat, MatStructure flag ) { ErrorCheck( ::KSPSetOperators( ksp, Amat, Pmat, flag ) ); }
213#else
214 inline static void KSPSetOperators (KSP ksp, Mat Amat, Mat Pmat ) { ErrorCheck( ::KSPSetOperators( ksp, Amat, Pmat ) ); }
215#endif
216 inline static void KSPSetTolerances ( KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits )
217 { ErrorCheck( ::KSPSetTolerances( ksp, rtol, abstol, dtol, maxits ) ); }
218 inline static void KSPSetInitialGuessNonzero( KSP ksp, PetscBool flg ) { ErrorCheck( ::KSPSetInitialGuessNonzero( ksp, flg ) ); };
219 inline static void KSPSetInitialGuessKnoll( KSP ksp, PetscBool flg ) { ErrorCheck( ::KSPSetInitialGuessKnoll( ksp, flg ) ); };
220 inline static void KSPSolve ( KSP ksp, Vec b, Vec x ) { ErrorCheck( ::KSPSolve( ksp, b, x ) ); }
221 inline static void KSPSetPC ( KSP ksp, PC pc ) { ErrorCheck( ::KSPSetPC( ksp, pc ) ); }
222
223 // preconditioning
224 template <class Comm>
225 inline static void PCCreate ( const Comm& comm, PC* pc)
226 {
227 ErrorCheck( ::PCCreate( castToPetscComm( comm ), pc ) );
228 }
229
230 inline static void PCDestroy ( PC* pc) {
231#if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
232 ErrorCheck( ::PCDestroy( *pc ) );
233#else
234 ErrorCheck( ::PCDestroy( pc ) );
235#endif
236 }
237 inline static void PCSetType ( PC pc, const PCType type ) { ErrorCheck( ::PCSetType( pc, type ) ); }
238 inline static void PCSetFromOptions ( PC pc ) { ErrorCheck( ::PCSetFromOptions( pc ) ); }
239 inline static void PCSetUp ( PC pc ) { ErrorCheck( ::PCSetUp( pc ) ); }
240 inline static void PCFactorSetLevels( PC pc, PetscInt level ) { ErrorCheck( ::PCFactorSetLevels( pc, level ) ); }
241 inline static void PCSORSetOmega( PC pc, PetscReal omega ) { ErrorCheck( ::PCSORSetOmega( pc, omega ) ); }
242 inline static void PCSORSetSymmetric( PC pc, MatSORType flag ) { ErrorCheck( ::PCSORSetSymmetric( pc, flag ) ); }
243#if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 9
244 inline static void PCFactorSetMatSolverPackage( PC pc, const MatSolverPackage type )
245 {
246 ErrorCheck( ::PCFactorSetMatSolverPackage( pc, type ) );
247 }
248#else
249 inline static void PCFactorSetMatSolverPackage( PC pc, const MatSolverType type )
250 {
251 ErrorCheck( ::PCFactorSetMatSolverType( pc, type ) );
252 }
253#endif
254 inline static void PCHYPRESetType( PC pc, const char* type )
255 {
256 ErrorCheck( ::PCHYPRESetType( pc, type ) );
257 }
258
259 // matrix routines
260 inline static void MatAssemblyBegin ( Mat mat, MatAssemblyType type ) { ErrorCheck( ::MatAssemblyBegin( mat, type ) ); }
261 inline static void MatAssemblyEnd ( Mat mat, MatAssemblyType type ) { ErrorCheck( ::MatAssemblyEnd( mat, type ) ); }
262 inline static void MatAssembled( Mat mat, PetscBool* assembled ) { ErrorCheck( ::MatAssembled( mat, assembled ) ); }
263
264 template <class Comm>
265 inline static void MatCreate ( const Comm& comm, Mat *A )
266 {
267 ErrorCheck( ::MatCreate( castToPetscComm( comm ), A) );
268 }
269
270 template <class Comm>
271 inline static void MatCreateBlockMat ( const Comm& comm,
272 Mat *A,
273 PetscInt m, PetscInt n, PetscInt bs, PetscInt nz, PetscInt* nnz )
274 {
275 ErrorCheck( ::MatCreateBlockMat( castToPetscComm( comm ), n, m, bs, nz, nnz, A) );
276 }
277 inline static void MatDestroy ( Mat *A ) {
278 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
279 ErrorCheck( ::MatDestroy( *A ) );
280 #else
281 ErrorCheck( ::MatDestroy( A ) );
282 #endif // PETSC_PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
283 }
284 inline static void MatSetUp( Mat mat )
285 {
286 ErrorCheck( ::MatSetUp(mat));
287 }
288 inline static void MatSetUp( Mat mat, PetscInt bs, PetscInt nz )
289 {
290 if (bs == 1)
291 {
292 ErrorCheck( ::MatSeqAIJSetPreallocation(mat,nz,PETSC_NULLPTR) );
293 ErrorCheck( ::MatMPIAIJSetPreallocation(mat,nz,PETSC_NULLPTR,nz/2,PETSC_NULLPTR) );
294 }
295 else
296 {
297 ErrorCheck( ::MatSeqBAIJSetPreallocation(mat,bs,nz,PETSC_NULLPTR) );
298 ErrorCheck( ::MatMPIBAIJSetPreallocation(mat,bs,nz,PETSC_NULLPTR,nz/2,PETSC_NULLPTR) );
299 }
300 // the following seems not to work for block matrix
301 ErrorCheck( ::MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) );
302 // the following only works for block matrix
303 // but should be used with MAT_NEW_NONZERO_LOCATIONS...
304 // ErrorCheck( ::MatSetOption(mat, MAT_USE_HASH_TABLE,PETSC_FALSE) );
305 ErrorCheck( ::MatSetUp(mat));
306 }
307 inline static void MatSetUp( Mat mat, PetscInt bs, const PetscInt *d_nnz, const PetscInt *o_nnz )
308 {
309 if (bs == 1)
310 {
311 ErrorCheck( ::MatSeqAIJSetPreallocation(mat,0,d_nnz ) );
312 ErrorCheck( ::MatMPIAIJSetPreallocation(mat,0,d_nnz,5,o_nnz) );
313 }
314 else
315 {
316 ErrorCheck( ::MatSeqBAIJSetPreallocation(mat,bs,0,d_nnz ) );
317 ErrorCheck( ::MatMPIBAIJSetPreallocation(mat,bs,0,d_nnz,5,PETSC_NULLPTR) );
318 }
319 // see previous comments
320 ErrorCheck( ::MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) );
321 ErrorCheck( ::MatSetUp(mat));
322 }
323
324 inline static void MatGetOwnershipRange ( Mat mat, PetscInt *m, PetscInt* n ) { ErrorCheck( ::MatGetOwnershipRange( mat, m, n ) ); }
325 inline static void MatGetSize ( Mat mat, PetscInt *m, PetscInt* n ) { ErrorCheck( ::MatGetSize( mat, m, n ) ); }
326 inline static void MatMult ( Mat mat, Vec x, Vec y ) { ErrorCheck( ::MatMult( mat, x, y ) ); }
327 inline static void MatSetBlockSize ( Mat A, PetscInt bs ) { ErrorCheck( ::MatSetBlockSize( A, bs ) ); }
328 inline static void MatSetSizes ( Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N ) { ErrorCheck( ::MatSetSizes( A, m, n, M, N ) ); }
329 inline static void MatSetFromOptions ( Mat B ) { ErrorCheck( ::MatSetFromOptions( B ) ); }
330 inline static void MatSetType ( Mat mat, const MatType matype ) { ErrorCheck( ::MatSetType( mat, matype ) ); }
331
332 inline static void MatSetValue ( Mat v, PetscInt i, PetscInt j, PetscScalar va, InsertMode mode )
333 {
334#ifdef _OPENMP
335#pragma omp critical
336#endif
337 {
338 ErrorCheck( ::MatSetValue( v, i, j, va, mode ) );
339 }
340 }
341
342 inline static void MatSetValues ( Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv )
343 {
344#ifdef _OPENMP
345#pragma omp critical
346#endif
347 {
348 ErrorCheck( ::MatSetValues( mat, m, idxm, n, idxn, v, addv ) );
349 }
350 }
351
352 inline static void MatSetValuesBlocked ( Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv )
353 {
354#ifdef _OPENMP
355#pragma omp critical
356#endif
357 {
358 ErrorCheck( ::MatSetValuesBlocked( mat, m, idxm, n, idxn, v, addv ) );
359 }
360 }
361
362 inline static void MatGetValues ( Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[] )
363 { ErrorCheck( ::MatGetValues( mat, m, idxm, n, idxn, v ) ); }
364
365 inline static void MatZeroRows ( Mat mat, PetscInt m, const PetscInt idxm[], const PetscScalar v )
366 {
367#ifdef _OPENMP
368#pragma omp critical
369#endif
370 {
371 ErrorCheck( ::MatZeroRows( mat, m, idxm, v, 0, 0 ) );
372 }
373 }
374
375 inline static void MatView ( Mat mat, PetscViewer viewer ) { ErrorCheck( ::MatView( mat, viewer ) ); }
376 inline static void MatZeroEntries ( Mat mat )
377 {
378#ifdef _OPENMP
379#pragma omp critical
380#endif
381 {
382 ErrorCheck( ::MatZeroEntries( mat ) );
383 }
384 }
385
386 inline static void PetscOptionsSetValue(const std::string& key, const std::string& value)
387 {
388 // add parameter to global data base (nullptr)
389 // TODO, add PetscOptions object to have several data bases
390 ErrorCheck( ::PetscOptionsSetValue(PETSC_NULLPTR, key.c_str(), value.c_str() ));
391 }
392
393 inline static void PetscOptionsClear()
394 {
395 // clear global data base (nullptr)
396 // TODO, add PetscOptions object to have several data bases
397 ErrorCheck( ::PetscOptionsClear(PETSC_NULLPTR) );
398 }
399
400 inline static void PetscBarrier ( PetscObject obj ) { ErrorCheck( ::PetscBarrier( obj ) ); }
401 inline static void PetscFinalize () { ErrorCheck( ::PetscFinalize() ); }
402 inline static void PetscInitialize( int *argc, char ***args, const char file[], const char help[] ) { ErrorCheck( ::PetscInitialize( argc, args, file, help ) ); }
403 inline static void PetscViewerASCIIOpen ( MPI_Comm comm, const char name[], PetscViewer *lab ) { ErrorCheck( ::PetscViewerASCIIOpen( comm, name, lab ) ); }
404 inline static void PetscViewerDestroy ( PetscViewer *viewer ) {
405 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
406 ErrorCheck( ::PetscViewerDestroy( *viewer ) );
407 #else
408 ErrorCheck( ::PetscViewerDestroy( viewer ) );
409 #endif
410 }
411 inline static void PetscViewerSetFormat ( PetscViewer viewer, PetscViewerFormat format )
412 {
413 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
414 ErrorCheck( ::PetscViewerSetFormat( viewer, format ) );
415 #else
416 ErrorCheck( ::PetscViewerPushFormat( viewer, format ) );
417 #endif
418 }
419 inline static void VecAssemblyBegin ( Vec vec ) { ErrorCheck( ::VecAssemblyBegin( vec ) ); }
420 inline static void VecAssemblyEnd ( Vec vec ) { ErrorCheck( ::VecAssemblyEnd( vec ) ); }
421 inline static void VecAXPY ( Vec y, PetscScalar alpha, Vec x) { ErrorCheck( ::VecAXPY( y, alpha, x ) ); }
422 inline static void VecCopy ( Vec x, Vec y ) { ErrorCheck( ::VecCopy( x, y ) ); }
423
424 template <class Comm>
425 inline static void VecCreate ( const Comm& comm, Vec *vec )
426 {
427 ErrorCheck( ::VecCreate( castToPetscComm( comm ), vec ) );
428 }
429
430 template <class Comm>
431 inline static void VecCreateGhost ( const Comm& comm, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[], Vec *vv )
432 { ErrorCheck( ::VecCreateGhost( castToPetscComm( comm ), n, N, nghost, ghosts, vv ) ); }
433
434 template <class Comm>
435 inline static void VecCreateGhostBlock ( const Comm& comm, PetscInt bs, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[], Vec *vv )
436 {
437#if PETSC_VERSION_MAJOR < 3 || (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR <= 2)
438 std::unique_ptr< PetscInt[] > ghostsCopy( new PetscInt[ nghost * bs ] );
439 for( PetscInt i = 0; i < nghost; ++i )
440 {
441 for( PetscInt j = 0; j < bs; ++j )
442 ghostsCopy[ i*bs + j ] = ghosts[ i ]*bs + j;
443 }
444 VecCreateGhost( n, N, nghost * bs, ghostsCopy.get(), vv );
445#else // #if PETSC_VERSION_MAJOR < 3 || (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR <= 2)
446 ErrorCheck( ::VecCreateGhostBlock( castToPetscComm( comm ), bs, n, N, nghost, ghosts, vv ) );
447#endif // #else // #if PETSC_VERSION_MAJOR < 3 || (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR <= 2)
448 }
449
450 inline static void VecDestroy ( Vec *v ) {
451 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
452 ErrorCheck( ::VecDestroy( *v ) );
453 #else
454 ErrorCheck( ::VecDestroy( v ) );
455 #endif // PETSC_PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 2
456 }
457 inline static void VecDot ( Vec x, Vec y, PetscScalar *val ) { ErrorCheck( ::VecDot( x, y, val ) ); }
458 inline static void VecDuplicate ( Vec v, Vec *newv ) { ErrorCheck( ::VecDuplicate( v, newv ) ); }
459 inline static void VecGetBlockSize ( Vec x, PetscInt* bs ) { ErrorCheck( ::VecGetBlockSize( x, bs ) ); }
460 inline static void VecGetLocalSize ( Vec x, PetscInt *size ) { ErrorCheck( ::VecGetLocalSize( x, size ) ); }
461 inline static void VecGetOwnershipRange ( Vec x, PetscInt *low, PetscInt *high ) { ErrorCheck( ::VecGetOwnershipRange( x, low, high ) ); }
462 inline static void VecGetSize ( Vec x, PetscInt *size ) { ErrorCheck( ::VecGetSize( x, size ) ); }
463 inline static void VecGetValues ( Vec x, PetscInt ni, const PetscInt ix[], PetscScalar y[] ) { ErrorCheck( ::VecGetValues( x, ni, ix, y ) ); }
464 inline static void VecGhostGetLocalForm ( Vec g, Vec *l ) { ErrorCheck( ::VecGhostGetLocalForm( g, l ) ); }
465 inline static void VecGhostRestoreLocalForm ( Vec g, Vec *l ) { ErrorCheck( ::VecGhostRestoreLocalForm( g, l ) ); }
466 inline static void VecGhostUpdateBegin ( Vec g, InsertMode insertmode, ScatterMode scattermode ) { ErrorCheck( ::VecGhostUpdateBegin( g, insertmode, scattermode ) ); }
467 inline static void VecGhostUpdateEnd ( Vec g, InsertMode insertmode, ScatterMode scattermode ) { ErrorCheck( ::VecGhostUpdateEnd( g, insertmode, scattermode ) ); }
468 inline static void VecNorm ( Vec x, NormType type, PetscReal *val ) { ErrorCheck( ::VecNorm( x, type, val ) ); }
469 inline static void VecScale ( Vec x, PetscScalar alpha ) { ErrorCheck( ::VecScale( x, alpha ) ); }
470 inline static void VecSet ( Vec x, PetscScalar alpha ) { ErrorCheck( ::VecSet( x, alpha ) ); }
471 inline static void VecSetBlockSize ( Vec x, PetscInt bs ) { ErrorCheck( ::VecSetBlockSize( x, bs ) ); }
472 inline static void VecSetFromOptions ( Vec vec ) { ErrorCheck( ::VecSetFromOptions( vec ) ); }
473 inline static void VecSetType ( Vec vec, const VecType method ) { ErrorCheck( ::VecSetType( vec, method ) ); }
474 inline static void VecSetSizes ( Vec v, PetscInt n, PetscInt N ) { ErrorCheck( ::VecSetSizes( v, n, N ) ); }
475 inline static void VecSetValue ( Vec v, int row, PetscScalar value, InsertMode mode ) { ErrorCheck( ::VecSetValue( v, row, value, mode ) ); }
476 inline static void VecSetValuesBlocked ( Vec v, PetscInt ni, const PetscInt xi[], const PetscScalar values[], InsertMode mode ) { ErrorCheck( ::VecSetValuesBlocked( v, ni, xi, values, mode ) ); }
477 inline static void VecView ( Vec vec, PetscViewer viewer ) { ErrorCheck( ::VecView( vec, viewer ) ); }
478
479 } // namespace Petsc
480
481} // namespace Dune
482
483#endif // #if HAVE_PETSC
484
485#endif // DUNE_FEM_PETSCCOMMON_HH
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
Implements an utility class that provides collective communication methods for sequential programs.
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr GeometryType line
GeometryType representing a line.
Definition: type.hh:498
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:117
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Standard Dune debug streams.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)