Actual source code: ex218.c


  2: static char help[] = "Tests MatShellTestMult()\n\n";

  4: #include <petscmat.h>

  6: typedef struct _n_User *User;
  7: struct _n_User {
  8:   Mat B;
  9: };

 11: static PetscErrorCode MatMult_User(Mat A, Vec X, Vec Y)
 12: {
 13:   User user;

 15:   PetscFunctionBegin;
 16:   PetscCall(MatShellGetContext(A, &user));
 17:   PetscCall(MatMult(user->B, X, Y));
 18:   PetscFunctionReturn(PETSC_SUCCESS);
 19: }

 21: static PetscErrorCode MatMultTranspose_User(Mat A, Vec X, Vec Y)
 22: {
 23:   User user;

 25:   PetscFunctionBegin;
 26:   PetscCall(MatShellGetContext(A, &user));
 27:   PetscCall(MatMultTranspose(user->B, X, Y));
 28:   PetscFunctionReturn(PETSC_SUCCESS);
 29: }

 31: static PetscErrorCode MyFunction(void *ctx, Vec x, Vec y)
 32: {
 33:   User user = (User)ctx;

 35:   PetscFunctionBegin;
 36:   PetscCall(MatMult(user->B, x, y));
 37:   PetscFunctionReturn(PETSC_SUCCESS);
 38: }

 40: int main(int argc, char **args)
 41: {
 42:   const PetscInt inds[]  = {0, 1};
 43:   PetscScalar    avals[] = {2, 3, 5, 7};
 44:   Mat            S;
 45:   User           user;
 46:   Vec            base;

 48:   PetscFunctionBeginUser;
 49:   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
 50:   PetscCall(PetscNew(&user));
 51:   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 2, 2, 2, NULL, &user->B));
 52:   PetscCall(MatSetUp(user->B));
 53:   PetscCall(MatSetValues(user->B, 2, inds, 2, inds, avals, INSERT_VALUES));
 54:   PetscCall(MatAssemblyBegin(user->B, MAT_FINAL_ASSEMBLY));
 55:   PetscCall(MatAssemblyEnd(user->B, MAT_FINAL_ASSEMBLY));
 56:   PetscCall(MatCreateVecs(user->B, &base, NULL));
 57:   PetscCall(MatCreateShell(PETSC_COMM_WORLD, 2, 2, 2, 2, user, &S));
 58:   PetscCall(MatSetUp(S));
 59:   PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MatMult_User));
 60:   PetscCall(MatShellSetOperation(S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_User));

 62:   PetscCall(MatShellTestMult(S, MyFunction, base, user, NULL));
 63:   PetscCall(MatShellTestMultTranspose(S, MyFunction, base, user, NULL));

 65:   PetscCall(VecDestroy(&base));
 66:   PetscCall(MatDestroy(&user->B));
 67:   PetscCall(MatDestroy(&S));
 68:   PetscCall(PetscFree(user));
 69:   PetscCall(PetscFinalize());
 70:   return 0;
 71: }

 73: /*TEST

 75:    test:
 76:      args: -mat_shell_test_mult_view -mat_shell_test_mult_transpose_view

 78: TEST*/