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