1 static char help[] = "Test MatMultHermitianTranspose() and MatMultHermitianTransposeAdd().\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc, char **args) { 6 Mat A, B, C; 7 Vec x, y, ys; 8 PetscInt i, j; 9 PetscScalar v; 10 PetscBool flg; 11 12 PetscFunctionBeginUser; 13 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 14 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 15 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 16 PetscCall(MatSetType(A, MATAIJ)); 17 PetscCall(MatSetFromOptions(A)); 18 PetscCall(MatSetUp(A)); 19 20 i = 0; 21 j = 0; 22 v = 2.0; 23 PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 24 i = 0; 25 j = 1; 26 v = 3.0 + 4.0 * PETSC_i; 27 PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 28 i = 1; 29 j = 0; 30 v = 5.0 + 6.0 * PETSC_i; 31 PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 32 i = 1; 33 j = 1; 34 v = 7.0 + 8.0 * PETSC_i; 35 PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 36 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 37 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 38 39 /* Create vectors */ 40 PetscCall(VecCreate(PETSC_COMM_WORLD, &y)); 41 PetscCall(VecSetSizes(y, PETSC_DECIDE, 2)); 42 PetscCall(VecSetFromOptions(y)); 43 PetscCall(VecDuplicate(y, &ys)); 44 PetscCall(VecDuplicate(y, &x)); 45 46 i = 0; 47 v = 10.0 + 11.0 * PETSC_i; 48 PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES)); 49 i = 1; 50 v = 100.0 + 120.0 * PETSC_i; 51 PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES)); 52 PetscCall(VecAssemblyBegin(x)); 53 PetscCall(VecAssemblyEnd(x)); 54 55 PetscCall(MatMultHermitianTranspose(A, x, y)); 56 PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 57 PetscCall(MatMultHermitianTransposeAdd(A, x, y, ys)); 58 PetscCall(VecView(ys, PETSC_VIEWER_STDOUT_WORLD)); 59 60 PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &B)); 61 PetscCall(MatCreateHermitianTranspose(A, &C)); 62 PetscCall(MatMultHermitianTransposeEqual(B, C, 4, &flg)); 63 PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "B^Hx != C^Hx"); 64 PetscCall(MatMultHermitianTransposeAddEqual(B, C, 4, &flg)); 65 PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "y+B^Hx != y+C^Hx"); 66 PetscCall(MatDestroy(&C)); 67 PetscCall(MatDestroy(&B)); 68 69 PetscCall(MatDestroy(&A)); 70 71 PetscCall(VecDestroy(&x)); 72 PetscCall(VecDestroy(&y)); 73 PetscCall(VecDestroy(&ys)); 74 PetscCall(PetscFinalize()); 75 return 0; 76 } 77 78 /*TEST 79 80 build: 81 requires: complex 82 test: 83 84 test: 85 suffix: 2 86 nsize: 2 87 88 TEST*/ 89