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