PetscWrapper.h 2.84 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
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
/******************************************************************************
 *
 * AMDiS - Adaptive multidimensional simulations
 *
 * Copyright (C) 2013 Dresden University of Technology. All Rights Reserved.
 * Web: https://fusionforge.zih.tu-dresden.de/projects/amdis
 *
 * Authors:
 * Simon Vey, Thomas Witkowski, Andreas Naumann, Simon Praetorius, et al.
 *
 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
 *
 *
 * This file is part of AMDiS
 *
 * See also license.opensource.txt in the distribution.
 *
 ******************************************************************************/



/** \file PetscWrapper.h */

#ifndef AMDIS_PETSCWRAPPER_H
#define AMDIS_PETSCWRAPPER_H

#include <mpi.h>
#include <petsc.h>

namespace AMDiS
{

  // functions for PETSc API changes
  namespace petsc
  {
    inline PetscErrorCode options_view(PetscViewer viewer)
    {
#if (PETSC_VERSION_MINOR >= 7)
      return PetscOptionsView(PETSC_NULL, viewer);
#else
      return PetscOptionsView(viewer);
#endif
    }

    inline PetscErrorCode options_insert_string(const char in_str[])
    {
#if (PETSC_VERSION_MINOR >= 7)
      return PetscOptionsInsertString(PETSC_NULL, in_str);
#else
      return PetscOptionsInsertString(in_str);
#endif
    }

    inline PetscErrorCode ksp_set_operators(KSP ksp, Mat Amat,Mat Pmat)
    {
#if (PETSC_VERSION_MINOR >= 5)
      return KSPSetOperators(ksp, Amat, Pmat);
#else
      return KSPSetOperators(ksp, Amat, Pmat, SAME_NONZERO_PATTERN);
#endif
    }

    inline PetscErrorCode ksp_get_operators(KSP ksp, Mat *Amat,Mat *Pmat)
    {
#if (PETSC_VERSION_MINOR >= 5)
      return KSPGetOperators(ksp, Amat, Pmat);
#else
      return KSPGetOperators(ksp, Amat, Pmat, SAME_NONZERO_PATTERN);
#endif
    }

    template <class Monitor>
    inline PetscErrorCode ksp_monitor_set(KSP ksp, Monitor monitor)
    {
#if (PETSC_VERSION_MINOR >= 7)
      PetscViewerAndFormat *vf;
      PetscErrorCode ierr;
      ierr = PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);CHKERRQ(ierr);
      ierr = KSPMonitorSet(ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
      return ierr;
#else
      return KSPMonitorSet(ksp, monitor, PETSC_NULL, PETSC_NULL);
#endif
    }

    inline PetscErrorCode mat_create_vecs(Mat mat,Vec *right,Vec *left)
    {
#if (PETSC_VERSION_MINOR >= 6)
      return MatCreateVecs(mat, right, left);
#else
      return MatGetVecs(mat, right, left);
#endif
    }

    inline PetscErrorCode mat_nullspace_remove(MatNullSpace sp,Vec vec)
    {
#if (PETSC_VERSION_MINOR >= 5)
      return MatNullSpaceRemove(sp, vec);
#else
      return MatNullSpaceRemove(sp, vec, PETSC_NULL);
#endif
    }

  } // end namespace petsc
} // end namespace AMDiS

#endif