mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Reduces CHKERRXX usage to where necessary.
The error checking macro makes it harder to read and harder to write, so instead we now only check for functions that can contain errors. Bounds and range checks are handled by PETSc and not OPM.
This commit is contained in:
@@ -128,22 +128,24 @@ namespace{
|
|||||||
PC preconditioner;
|
PC preconditioner;
|
||||||
|
|
||||||
OEM_DATA( const int size ) {
|
OEM_DATA( const int size ) {
|
||||||
CHKERRXX( VecCreate( PETSC_COMM_WORLD, &solution ) );
|
VecCreate( PETSC_COMM_WORLD, &solution );
|
||||||
CHKERRXX( VecSetSizes( solution, PETSC_DECIDE, size ) );
|
auto err = VecSetSizes( solution, PETSC_DECIDE, size );
|
||||||
CHKERRXX( VecSetFromOptions( solution ) );
|
CHKERRXX( err );
|
||||||
CHKERRXX( VecDuplicate( solution, &rhs ) );
|
VecSetFromOptions( solution );
|
||||||
|
VecDuplicate( solution, &rhs );
|
||||||
|
|
||||||
CHKERRXX( MatCreate( PETSC_COMM_WORLD, &A ) );
|
MatCreate( PETSC_COMM_WORLD, &A );
|
||||||
CHKERRXX( MatSetSizes( A, PETSC_DECIDE, PETSC_DECIDE, size, size ) );
|
err = MatSetSizes( A, PETSC_DECIDE, PETSC_DECIDE, size, size );
|
||||||
CHKERRXX( MatSetFromOptions( A ) );
|
CHKERRXX( err );
|
||||||
CHKERRXX( MatSetUp( A ) );
|
MatSetFromOptions( A );
|
||||||
|
MatSetUp( A );
|
||||||
}
|
}
|
||||||
|
|
||||||
~OEM_DATA() {
|
~OEM_DATA() {
|
||||||
CHKERRXX( VecDestroy( &rhs ) );
|
VecDestroy( &rhs );
|
||||||
CHKERRXX( VecDestroy( &solution ) );
|
VecDestroy( &solution );
|
||||||
CHKERRXX( MatDestroy( &A ) );
|
MatDestroy( &A );
|
||||||
CHKERRXX( KSPDestroy( &ksp ) );
|
KSPDestroy( &ksp );
|
||||||
}
|
}
|
||||||
};
|
};
|
||||||
|
|
||||||
@@ -154,11 +156,11 @@ namespace{
|
|||||||
PetscScalar* vec;
|
PetscScalar* vec;
|
||||||
PetscInt size;
|
PetscInt size;
|
||||||
|
|
||||||
CHKERRXX( VecGetLocalSize( v, &size ) );
|
VecGetLocalSize( v, &size );
|
||||||
CHKERRXX( VecGetArray( v, &vec ) );
|
VecGetArray( v, &vec );
|
||||||
|
|
||||||
std::memcpy( vec, x, size * sizeof( double ) );
|
std::memcpy( vec, x, size * sizeof( double ) );
|
||||||
CHKERRXX( VecRestoreArray( v, &vec ) );
|
VecRestoreArray( v, &vec );
|
||||||
}
|
}
|
||||||
|
|
||||||
void from_petsc_vec( double* x, Vec v ) {
|
void from_petsc_vec( double* x, Vec v ) {
|
||||||
@@ -168,11 +170,11 @@ namespace{
|
|||||||
PetscScalar* vec;
|
PetscScalar* vec;
|
||||||
PetscInt size;
|
PetscInt size;
|
||||||
|
|
||||||
CHKERRXX( VecGetLocalSize( v, &size ) );
|
VecGetLocalSize( v, &size );
|
||||||
CHKERRXX( VecGetArray( v, &vec ) );
|
VecGetArray( v, &vec );
|
||||||
|
|
||||||
std::memcpy( x, vec, size * sizeof( double ) );
|
std::memcpy( x, vec, size * sizeof( double ) );
|
||||||
CHKERRXX( VecRestoreArray( v, &vec ) );
|
VecRestoreArray( v, &vec );
|
||||||
}
|
}
|
||||||
|
|
||||||
void to_petsc_mat( const int size, const int nonzeros,
|
void to_petsc_mat( const int size, const int nonzeros,
|
||||||
@@ -183,11 +185,14 @@ namespace{
|
|||||||
|
|
||||||
if( nrows == 0 ) continue;
|
if( nrows == 0 ) continue;
|
||||||
|
|
||||||
for( int j = ia[ i ]; j < ia[ i + 1 ]; ++j )
|
for( int j = ia[ i ]; j < ia[ i + 1 ]; ++j ) {
|
||||||
CHKERRXX( MatSetValues( A, 1, &i, 1, &ja[ j ], &sa[ j ], INSERT_VALUES ) );
|
auto err = MatSetValues( A, 1, &i, 1, &ja[ j ], &sa[ j ], INSERT_VALUES );
|
||||||
|
CHKERRXX( err );
|
||||||
|
}
|
||||||
|
|
||||||
CHKERRXX( MatAssemblyBegin( A, MAT_FINAL_ASSEMBLY ) );
|
auto err = MatAssemblyBegin( A, MAT_FINAL_ASSEMBLY );
|
||||||
CHKERRXX( MatAssemblyEnd( A, MAT_FINAL_ASSEMBLY ) );
|
CHKERRXX( err );
|
||||||
|
MatAssemblyEnd( A, MAT_FINAL_ASSEMBLY );
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
@@ -198,23 +203,28 @@ namespace{
|
|||||||
PetscReal residual;
|
PetscReal residual;
|
||||||
KSPConvergedReason reason;
|
KSPConvergedReason reason;
|
||||||
|
|
||||||
CHKERRXX( KSPCreate( PETSC_COMM_WORLD, &t.ksp ) );
|
KSPCreate( PETSC_COMM_WORLD, &t.ksp );
|
||||||
CHKERRXX( KSPSetOperators( t.ksp, t.A, t.A, DIFFERENT_NONZERO_PATTERN ) );
|
KSPSetOperators( t.ksp, t.A, t.A, DIFFERENT_NONZERO_PATTERN );
|
||||||
CHKERRXX( KSPGetPC( t.ksp, &t.preconditioner ) );
|
KSPGetPC( t.ksp, &t.preconditioner );
|
||||||
CHKERRXX( KSPSetType( t.ksp, method ) );
|
auto err = KSPSetType( t.ksp, method );
|
||||||
CHKERRXX( PCSetType( t.preconditioner, pcname ) );
|
CHKERRXX( err );
|
||||||
CHKERRXX( KSPSetTolerances( t.ksp, rtol, atol, dtol, maxits ) );
|
err = PCSetType( t.preconditioner, pcname );
|
||||||
CHKERRXX( KSPSetFromOptions( t.ksp ) );
|
CHKERRXX( err );
|
||||||
CHKERRXX( KSPSetInitialGuessNonzero( t.ksp, PETSC_TRUE ) );
|
err = KSPSetTolerances( t.ksp, rtol, atol, dtol, maxits );
|
||||||
CHKERRXX( KSPSolve( t.ksp, t.rhs, t.solution ) );
|
CHKERRXX( err );
|
||||||
CHKERRXX( KSPGetConvergedReason( t.ksp, &reason ) );
|
err = KSPSetFromOptions( t.ksp );
|
||||||
CHKERRXX( KSPGetIterationNumber( t.ksp, &its ) );
|
CHKERRXX( err );
|
||||||
CHKERRXX( KSPGetResidualNorm( t.ksp, &residual ) );
|
KSPSetInitialGuessNonzero( t.ksp, PETSC_TRUE );
|
||||||
|
KSPSolve( t.ksp, t.rhs, t.solution );
|
||||||
|
KSPGetConvergedReason( t.ksp, &reason );
|
||||||
|
KSPGetIterationNumber( t.ksp, &its );
|
||||||
|
KSPGetResidualNorm( t.ksp, &residual );
|
||||||
|
|
||||||
if( ksp_view )
|
if( ksp_view )
|
||||||
CHKERRXX( KSPView( t.ksp, PETSC_VIEWER_STDOUT_WORLD ) );
|
KSPView( t.ksp, PETSC_VIEWER_STDOUT_WORLD );
|
||||||
|
|
||||||
CHKERRXX( PetscPrintf( PETSC_COMM_WORLD, "KSP Iterations %D, Final Residual %G\n", its, residual ) );
|
err = PetscPrintf( PETSC_COMM_WORLD, "KSP Iterations %D, Final Residual %G\n", its, residual );
|
||||||
|
CHKERRXX( err );
|
||||||
}
|
}
|
||||||
|
|
||||||
} // anonymous namespace.
|
} // anonymous namespace.
|
||||||
|
|||||||
Reference in New Issue
Block a user