[petsc-users] problem with MatShellGetContext

Matthew Knepley knepley at gmail.com
Wed Aug 5 06:38:20 CDT 2015


On Wed, Aug 5, 2015 at 4:15 AM, Nicolas Pozin <nicolas.pozin at inria.fr>
wrote:

> Hello,
>
> I'm trying to solve a system with a matrix free operator and through
> conjugate gradient method.
> To make ideas clear, I set up the following simple example (I am using
> petsc-3.6) and I get this error message :
>

Yes, you are passing a C++ function userMult, so the compiler sticks "this"
in as the first argument. We do not
recommend this kind of wrapping.

  Thanks,

    Matt


> "
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Invalid argument!
> [0]PETSC ERROR: Wrong type of object: Parameter # 1!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.4.3, Oct, 15, 2013
> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [0]PETSC ERROR: See docs/index.html for manual pages.
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: ./test on a ubuntu_release named pl-59080 by npozin Wed
> Aug  5 10:55:26 2015
> [0]PETSC ERROR: Libraries linked from
> /home/npozin/Felisce_libraries/petsc_3.4.3/ubuntu_release/lib
> [0]PETSC ERROR: Configure run at Wed Jul 22 16:18:36 2015
> [0]PETSC ERROR: Configure options PETSC_ARCH=ubuntu_release --with-cxx=g++
> --with-fc=gfortran --with-cc=gcc --with-x=0 --download-openmpi
> --download-f-blas-lapack --download-superlu --download-superlu_dist
> --with-superlu_dist=1 --download-metis --download-mumps --download-parmetis
> --with-superlu_dist=1 --download-boost --with-boost=1 --download-scalapack
> with-external-packages-dir=/home/npozin/Felisce_libraries/petsc_3.4.3/packages
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: MatShellGetContext() line 202 in
> /home/npozin/Felisce_libraries/petsc_3.4.3/src/mat/impls/shell/shell.c
> End userMult
> [0]PETSC ERROR: MatMult() line 2179 in
> /home/npozin/Felisce_libraries/petsc_3.4.3/src/mat/interface/matrix.c
> [0]PETSC ERROR: KSP_MatMult() line 204 in
> /home/npozin/Felisce_libraries/petsc_3.4.3/include/petsc-private/kspimpl.h
> [0]PETSC ERROR: KSPSolve_CG() line 219 in
> /home/npozin/Felisce_libraries/petsc_3.4.3/src/ksp/ksp/impls/cg/cg.c
> [0]PETSC ERROR: KSPSolve() line 441 in
> /home/npozin/Felisce_libraries/petsc_3.4.3/src/ksp/ksp/interface/itfunc.c
> "
>
> I don't understand where the problem comes from with the matrix argument
> of MatShellGetContext.
> Any idea on what I do wrong?
>
> Thanks a lot,
> Nicolas
>
>
>
> #include <iostream>
> #include <petscksp.h>
>
> using namespace std;
>
>
> typedef struct {
>   int val;
> } MyCtx;
>
>
> class ShellClass {
>   Mat matShell;
>   KSP ksp;
>   PC pc;
>   Vec x;
>   Vec b;
>
> public:
>   void userMult(Mat Amat, Vec x, Vec y) {
>     cout << "Inside userMult" << endl;
>
>     MyCtx *ctx;
>     MatShellGetContext(Amat, (void *) ctx);
>
>     cout << "End userMult" << endl;
>   }
>
>   void solveShell() {
>     // context
>     MyCtx *ctx = new MyCtx;
>     ctx->val = 42;
>
>     // pc
>     PCCreate(PETSC_COMM_WORLD, &pc);
>     PCSetType(pc, PCNONE);
>
>     // ksp
>     KSPCreate(PETSC_COMM_WORLD, &ksp);
>     KSPSetType(ksp, KSPCG);
>     KSPSetPC(ksp, pc);
>     KSPSetFromOptions(ksp);
>
>     // matshell
>     int m = 10;
>     int n = 10;
>     MatCreateShell(PETSC_COMM_WORLD, m, n, PETSC_DETERMINE,
> PETSC_DETERMINE, ctx, &matShell);
>     MatShellSetOperation(matShell, MATOP_MULT,
> (void(*)(void))&ShellClass::userMult);
>
>
>     // create vectors
>     MatCreateVecs(matShell, &x, 0);
>     VecDuplicate(x, &b);
>     VecSet(b, 1.);
>
>     // set operators
>     KSPSetOperators(ksp, matShell, matShell);
>
>     // solve (call to userMult)
>     KSPSolve(ksp, b, x);
>   }
> };
>
>
>
> int main(int argc, char** argv) {
>   PetscInitialize(&argc, &argv, NULL, NULL);
>
>   ShellClass foo;
>   foo.solveShell();
>
>   PetscFinalize();
>   return 0;
> }
>
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150805/c947aec7/attachment.html>


More information about the petsc-users mailing list