xref: /petsc/src/mat/tests/ex175.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateHermitianTranspose().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat         C,C_htransposed,Cht,C_empty;
9c4762a1bSJed Brown   PetscInt    i,j,m = 10,n = 10;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown 
12*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
13c4762a1bSJed Brown   /* Create a complex non-hermitian matrix */
145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
18c4762a1bSJed Brown   for (i=0; i<m; i++) {
19c4762a1bSJed Brown     for (j=0; j<n; j++) {
20c4762a1bSJed Brown       v = 0.0 - 1.0*PETSC_i;
215f80ce2aSJacob Faibussowitsch       if (i>j && i-j<2)   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
22c4762a1bSJed Brown     }
23c4762a1bSJed Brown   }
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
26c4762a1bSJed Brown 
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHermitianTranspose(C, &C_htransposed));
28c4762a1bSJed Brown 
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_SELF));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C_htransposed,MAT_COPY_VALUES,&Cht));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(Cht,PETSC_VIEWER_STDOUT_SELF));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C_htransposed,MAT_DO_NOT_COPY_VALUES,&C_empty));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C_empty,PETSC_VIEWER_STDOUT_SELF));
34c4762a1bSJed Brown 
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C_htransposed));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Cht));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C_empty));
39*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
40*b122ec5aSJacob Faibussowitsch   return 0;
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown /*TEST
44c4762a1bSJed Brown 
45c4762a1bSJed Brown    build:
46c4762a1bSJed Brown      requires: complex
47c4762a1bSJed Brown 
48c4762a1bSJed Brown    test:
49c4762a1bSJed Brown      output_file: output/ex175.out
50c4762a1bSJed Brown 
51c4762a1bSJed Brown TEST*/
52