1 2 static char help[] = "Tests MatCreateHermitianTranspose().\n\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat C,C_htransposed,Cht,C_empty; 9 PetscInt i,j,m = 10,n = 10; 10 PetscScalar v; 11 12 PetscFunctionBeginUser; 13 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 14 /* Create a complex non-hermitian matrix */ 15 PetscCall(MatCreate(PETSC_COMM_SELF,&C)); 16 PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n)); 17 PetscCall(MatSetFromOptions(C)); 18 PetscCall(MatSetUp(C)); 19 for (i=0; i<m; i++) { 20 for (j=0; j<n; j++) { 21 v = 0.0 - 1.0*PETSC_i; 22 if (i>j && i-j<2) PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 23 } 24 } 25 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 26 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 27 28 PetscCall(MatCreateHermitianTranspose(C, &C_htransposed)); 29 30 PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF)); 31 PetscCall(MatDuplicate(C_htransposed,MAT_COPY_VALUES,&Cht)); 32 PetscCall(MatView(Cht,PETSC_VIEWER_STDOUT_SELF)); 33 PetscCall(MatDuplicate(C_htransposed,MAT_DO_NOT_COPY_VALUES,&C_empty)); 34 PetscCall(MatView(C_empty,PETSC_VIEWER_STDOUT_SELF)); 35 36 PetscCall(MatDestroy(&C)); 37 PetscCall(MatDestroy(&C_htransposed)); 38 PetscCall(MatDestroy(&Cht)); 39 PetscCall(MatDestroy(&C_empty)); 40 PetscCall(PetscFinalize()); 41 return 0; 42 } 43 44 /*TEST 45 46 build: 47 requires: complex 48 49 test: 50 output_file: output/ex175.out 51 52 TEST*/ 53