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 PetscErrorCode ierr; 10 PetscInt i,j; 11 PetscScalar v; 12 PetscBool flg; 13 14 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 15 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 16 ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 17 ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 18 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 19 ierr = MatSetUp(A);CHKERRQ(ierr); 20 21 i = 0; j = 0; v = 2.0; 22 ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 23 i = 0; j = 1; v = 3.0 + 4.0*PETSC_i; 24 ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 25 i = 1; j = 0; v = 5.0 + 6.0*PETSC_i; 26 ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 27 i = 1; j = 1; v = 7.0 + 8.0*PETSC_i; 28 ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 29 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 30 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 31 32 /* Create vectors */ 33 ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr); 34 ierr = VecSetSizes(y,PETSC_DECIDE,2);CHKERRQ(ierr); 35 ierr = VecSetFromOptions(y);CHKERRQ(ierr); 36 ierr = VecDuplicate(y,&ys);CHKERRQ(ierr); 37 ierr = VecDuplicate(y,&x);CHKERRQ(ierr); 38 39 i = 0; v = 10.0 + 11.0*PETSC_i; 40 ierr = VecSetValues(x,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 41 i = 1; v = 100.0 + 120.0*PETSC_i; 42 ierr = VecSetValues(x,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 43 ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 44 ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 45 46 ierr = MatMultHermitianTranspose(A,x,y);CHKERRQ(ierr); 47 ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 48 ierr = MatMultHermitianTransposeAdd(A,x,y,ys);CHKERRQ(ierr); 49 ierr = VecView(ys,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 50 51 ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 52 ierr = MatCreateHermitianTranspose(A,&C);CHKERRQ(ierr); 53 ierr = MatMultHermitianTransposeEqual(B,C,4,&flg);CHKERRQ(ierr); 54 PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"B^Hx != C^Hx"); 55 ierr = MatMultHermitianTransposeAddEqual(B,C,4,&flg);CHKERRQ(ierr); 56 PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+B^Hx != y+C^Hx"); 57 ierr = MatDestroy(&C);CHKERRQ(ierr); 58 ierr = MatDestroy(&B);CHKERRQ(ierr); 59 60 ierr = MatDestroy(&A);CHKERRQ(ierr); 61 62 ierr = VecDestroy(&x);CHKERRQ(ierr); 63 ierr = VecDestroy(&y);CHKERRQ(ierr); 64 ierr = VecDestroy(&ys);CHKERRQ(ierr); 65 ierr = PetscFinalize(); 66 return ierr; 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