Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Stefano's fixes and improvements #2

Draft
wants to merge 3 commits into
base: master
Choose a base branch
from
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions .gitignore
Original file line number Diff line number Diff line change
@@ -0,0 +1 @@
build/
8 changes: 0 additions & 8 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,6 @@ project(ngspetsc)

cmake_minimum_required(VERSION 3.8)

option (PETSC_WITH_HYPRE "Does PETSc have hypre enabled?" ON)
option (PETSC_COMPLEX "Are we building the complex version of the interface?" OFF)

set(CMAKE_MODULE_PATH "${CMAKE_MODULE_PATH}" "${CMAKE_CURRENT_SOURCE_DIR}/cmake/cmake_modules")
Expand Down Expand Up @@ -32,12 +31,5 @@ else (NOT PETSC_COMPLEX)
set (MODULE_NAME "ngs_petsc_complex")
endif (NOT PETSC_COMPLEX)

if(NOT PETSC_WITH_HYPRE)
message(STATUS "Turning OFF Hypre Auxiliary space Preconditioners.")
set(CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} -DPETSC_NO_HYPRE")
else(NOT PETSC_WITH_HYPRE)
message(STATUS "Turning ON Hypre Auxiliary space Preconditioners.")
endif(NOT PETSC_WITH_HYPRE)

add_subdirectory(src)
add_subdirectory(python)
5 changes: 0 additions & 5 deletions README.rst
Original file line number Diff line number Diff line change
Expand Up @@ -11,8 +11,3 @@ override this with

cmake -DPETSC_EXECUTABLE_RUNS=YES

If your PETSc installation has beed configured without hypre, you have to tell cmake:


cmake -DPETSC_WITH_HYPRE=OFF

5 changes: 5 additions & 0 deletions src/petsc_interface.hpp
Original file line number Diff line number Diff line change
@@ -1,6 +1,11 @@

#include <comp.hpp>
// SZ
// consider removing this include from here
#include "petsc.h"

// SZ
// why you need python stuff in a C++ interface declaration?
#include <python_ngstd.hpp>

#include "typedefs.hpp"
Expand Down
4 changes: 0 additions & 4 deletions src/petsc_ksp.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -98,12 +98,8 @@ namespace ngs_petsc_interface
KSPSetUp(GetKSP());

petsc_rhs = GetMatrix()->GetRowMap()->CreatePETScVector();
// VecAssemblyBegin(petsc_rhs);
// VecAssemblyEnd(petsc_rhs);

petsc_sol = GetMatrix()->GetColMap()->CreatePETScVector();
// VecAssemblyBegin(petsc_sol);
// VecAssemblyEnd(petsc_sol);
}


Expand Down
5 changes: 5 additions & 0 deletions src/petsc_ksp.hpp
Original file line number Diff line number Diff line change
@@ -1,6 +1,11 @@
#ifndef FILE_NGSPETSC_KSP_HPP
#define FILE_NGSPETSC_KSP_HPP

// SZ
// Each specific hpp file should be able to resolve all the symbols defined in it
// So, at least for PETSc stuff, you can either include typedefs.hpp or forward declare the needed symbols (KSP, Mat, PC, etc)
// Since this is a possibly public header, you should avoid including petscksp.h

namespace ngs_petsc_interface
{

Expand Down
12 changes: 12 additions & 0 deletions src/petsc_linalg.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -125,6 +125,7 @@ namespace ngs_petsc_interface
template<class TM>
PETScMat CreatePETScMatSeqBAIJFromSymmetric (shared_ptr<ngs::SparseMatrixSymmetric<TM>> spmat, shared_ptr<ngs::BitArray> rss, shared_ptr<ngs::BitArray> css)
{
// SZ: it seems you could interface to SBAIJ too
static ngs::Timer t(string("CreatePETScMatSeqBAIJFromSymmetric<Mat<") + to_string(ngs::mat_traits<TM>::HEIGHT) + string(">>")); ngs::RegionTimer rt(t);

// row map (map for a row)
Expand Down Expand Up @@ -157,6 +158,8 @@ namespace ngs_petsc_interface
}
PETScMat petsc_mat;
if (bh == 1)
// SZ: here and in all other similar calls: if you now the nonzero per row, then pass 0 (not -1) as
// total number of nonzeros
{ MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, -1, &nzepr[1], &petsc_mat); }
else
{ MatCreateSeqBAIJ(PETSC_COMM_SELF, bh, nrows, ncols, -1, &nzepr[1], &petsc_mat); }
Expand Down Expand Up @@ -299,6 +302,13 @@ namespace ngs_petsc_interface
} // CreatePETScMatSeqBAIJ


// SZ
// I would not use BAIJ by default if the block size > 2
// AIJ is much more supported, BAIJ does not really deliver impressive performances improvements
// Some solvers (GAMG, BDDC, FIELDSPLIT) use block size information but this information
// is valid for AIJ too, if you use MatSetBlockSize(aijmat,bs);
// My suggestion here is to always default to AIJ and provide an extra argument to this
// function to instead use blocked formats. BTW, It seems SBAIJ support is missing, am I correct?
PETScMat CreatePETScMatSeq (shared_ptr<ngs::BaseMatrix> mat, shared_ptr<ngs::BitArray> rss, shared_ptr<ngs::BitArray> css)
{
if (auto spm = dynamic_pointer_cast<ngs::SparseMatrixTM<PETScScalar>>(mat))
Expand Down Expand Up @@ -758,6 +768,8 @@ namespace ngs_petsc_interface
static ngs::Timer t("PETScMatrix::UpdateValues"); ngs::RegionTimer rt(t);

// TODO: if we have converted the matrix from BAIJ to AIJ, is SetValuesBlocked inefficient??
// SZ: answer: it shouldn't be inefficient: it expands by blocksize
// the blocked set of indices and then call MatSetValues
auto parmat = dynamic_pointer_cast<ngs::ParallelMatrix>(ngs_mat);
shared_ptr<BaseMatrix> mat = (parmat == nullptr) ? ngs_mat : parmat->GetMatrix();

Expand Down
3 changes: 3 additions & 0 deletions src/petsc_linalg.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,8 @@
/**
Matrices / Vectors
**/
// SZ
// See comments about header inclusions in petsc_ksp.hpp

namespace ngs_petsc_interface
{
Expand Down Expand Up @@ -142,6 +144,7 @@ namespace ngs_petsc_interface

Can be any kind of BaseMatrix
**/
// SZ An alternative is to have another Enum value in MAT_TYPE, namely SHELL, and you do not need a special class for it
class FlatPETScMatrix : public PETScBaseMatrix
{
public:
Expand Down
7 changes: 6 additions & 1 deletion src/petsc_ops.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -130,17 +130,22 @@ namespace ngs_petsc_interface
void Ngs2PETSc (shared_ptr<ngs::BaseVector> ngs_vec, ::Vec petsc_vec)
{
ngs_vec->Cumulate();
VecAssemblyBegin(petsc_vec);
auto fv = ngs_vec->FVDouble();
for (auto k : Range(loc))
buf[k] = fv(loc_inds[k]);
VecSetValues(petsc_vec, loc, &glob_nums[0], &buf[0], INSERT_VALUES);
VecAssemblyBegin(petsc_vec);
VecAssemblyEnd(petsc_vec);
}
// petsc -> ngsolve
void PETSc2Ngs (shared_ptr<ngs::BaseVector> ngs_vec, ::Vec petsc_vec)
{
ngs_vec->Distribute();
// SZ
// This is odd, and rarely used from PETSc. Also, off-processor retrieval is disabled
// Can you use VecGetArrayRead/VecRestoreArrayRead semantics?
// In alternative, if Ngs2PETSc and PETSc2NGs needs off-process data, you can use
// the VecScatter object
VecGetValues(petsc_vec, loc, &glob_nums[0], &buf[0]);
auto fv = ngs_vec->FVDouble();
fv = 0.0;
Expand Down
8 changes: 5 additions & 3 deletions src/petsc_pc.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -189,7 +189,7 @@ namespace ngs_petsc_interface
}


#ifdef PETSC_HAS_HYPRE
#ifdef PETSC_HAVE_HYPRE
PETScHypreAuxiliarySpacePC :: PETScHypreAuxiliarySpacePC (shared_ptr<ngs::BilinearForm> _bfa, const ngs::Flags & _aflags, const string _aname)
: PETSc2NGsPrecond(_bfa, _aflags, _aname)
{
Expand Down Expand Up @@ -455,6 +455,8 @@ namespace ngs_petsc_interface

// Now set the PCs we already had before
KSP* ksps; PetscInt n;
// SZ: be careful with this function, since those KSPs are different depending on the
// fieldsplit type, see the online documentation
PCFieldSplitGetSubKSP(GetPETScPC(), &n, &ksps);
for (auto k : Range(fields.Size())) {
auto field = fields[k];
Expand All @@ -471,7 +473,7 @@ namespace ngs_petsc_interface

ngs::RegisterPreconditioner<PETSc2NGsPrecond> registerPETSc2NGsPrecond("petsc_pc");

#ifdef PETSC_HAS_HYPRE
#ifdef PETSC_HAVE_HYPRE
ngs::RegisterPreconditioner<PETScHypreAMS> registerPETScHypreAMS("petsc_pc_hypre_ams");
#endif

Expand Down Expand Up @@ -513,7 +515,7 @@ namespace ngs_petsc_interface
.def("Finalize", [](shared_ptr<PETSc2NGsPrecond> & pc)
{ pc->FinalizeLevel(); } );

#ifdef PETSC_HAS_HYPRE
#ifdef PETSC_HAVE_HYPRE
py::class_<PETScHypreAuxiliarySpacePC, shared_ptr<PETScHypreAuxiliarySpacePC>, PETSc2NGsPrecond>
(m, "HypreAuxiliarySpacePrecond", "ADS/AMS from hypre package")
.def ("SetGradientMatrix" , [](shared_ptr<PETScHypreAuxiliarySpacePC> & pc, shared_ptr<PETScMatrix> mat)
Expand Down
4 changes: 3 additions & 1 deletion src/petsc_pc.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,8 @@ namespace ngs_petsc_interface
// };


// SZ
// Why do you inherit from FlatPETScMatrix?
/** An NGSolve-BaseMatrix, wrapped to PETSc as a PC **/
class NGs2PETScPrecond : public PETScBasePrecond,
public FlatPETScMatrix
Expand Down Expand Up @@ -104,7 +106,7 @@ namespace ngs_petsc_interface
Array<shared_ptr<NGs2PETScPrecond>> keep_alive;
};

#ifdef PETSC_HAS_HYPRE
#ifdef PETSC_HAVE_HYPRE

class PETScHypreAuxiliarySpacePC : public PETSc2NGsPrecond
{
Expand Down
12 changes: 11 additions & 1 deletion src/petsc_snes.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -101,6 +101,7 @@ namespace ngs_petsc_interface
{ SetOptions (_opts, "", NULL); }

// Create Vector to hold F(x)
// SZ: this not strictly needed, since Petsc will create a Vec for you if you pass NULL
func_vec = jac_mat->GetRowMap()->CreatePETScVector();

// Set (non-linear) function evaluation, f = F(x)
Expand Down Expand Up @@ -234,9 +235,18 @@ namespace ngs_petsc_interface
{
auto& self = *( (PETScSNES*) ctx);
HeapReset hr(*self.use_lh);


// SZ: (do you really want users reading WTF in your code???)
// actually, this happens if snes_mf or snes_mf_operator are set
// if(A != self.jac_mat->GetPETScMat()) { throw Exception("A != JAC BUFFER, WTF?"); }

// SZ: B is the correct choice here. However, I think to cover all possible matrix-free
// options, you need to add the following 4 lines (this will reset the shell matrix to its proper state)
//
//if (A && A != B) {
// ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
// ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
//}
if(B != self.jac_mat->GetPETScMat()) { throw Exception("B != JAC BUFFER, WTF?"); }

self.jac_mat->GetRowMap()->PETSc2NGs(*self.lin_vec, x);
Expand Down
6 changes: 2 additions & 4 deletions src/typedefs.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -4,10 +4,6 @@
namespace ngs_petsc_interface
{

#ifndef PETSC_NO_HYPRE
#define PETSC_HAS_HYPRE
#endif

namespace ngs = ngcomp;
using ngs::Array;
using ngs::Range;
Expand All @@ -32,6 +28,8 @@ namespace ngs_petsc_interface
#else
static_assert( is_same<PetscScalar, double>::value, "Trying to compile the real interface with a complex PETSc installation, (set -DPETSC_COMPLEX=OFF)!");
#endif
//SZ
// Not sure if you support PETSc compiled with 64bit indices

// static_assert( (is_same<PetscScalar, double>::value || is_same<PetscScalar, ngs::Complex>::value), "Need double or complex PETSc version!");
// static_assert( (is_same<PetscScalar, double>::value), "Not a double PETSc version!");
Expand Down