Actual source code: ex175.c
petsc-3.7.3 2016-08-01
2: static char help[] = "Tests MatCreateHermitianTranspose().\n\n";
4: #include <petscmat.h>
8: int main(int argc,char **args)
9: {
10: Mat C,C_htransposed,Cht,C_empty;
11: PetscInt i,j,m = 10,n = 10;
13: PetscScalar v;
15: PetscInitialize(&argc,&args,(char*)0,help);
17: /* Create a complex non-hermitian matrix */
18: MatCreate(PETSC_COMM_SELF,&C);
19: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n);
20: MatSetFromOptions(C);
21: MatSetUp(C);
22: for (i=0; i<m; i++) {
23: for (j=0; j<n; j++) {
24: v = 0.0 - 1.0*PETSC_i;
25: if (i>j && i-j<2) {MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);}
26: }
27: }
28: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
29: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
31: MatCreateHermitianTranspose(C, &C_htransposed);
33: MatView(C,PETSC_VIEWER_STDOUT_SELF);
34: MatDuplicate(C_htransposed,MAT_COPY_VALUES,&Cht);
35: MatView(Cht,PETSC_VIEWER_STDOUT_SELF);
36: MatDuplicate(C_htransposed,MAT_DO_NOT_COPY_VALUES,&C_empty);
37: MatView(C_empty,PETSC_VIEWER_STDOUT_SELF);
39: MatDestroy(&C);
40: MatDestroy(&C_htransposed);
41: MatDestroy(&Cht);
42: MatDestroy(&C_empty);
43: PetscFinalize();
44: return 0;
45: }