xref: /petsc/src/mat/tests/ex175.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
14c4762a1bSJed Brown   /* Create a complex non-hermitian matrix */
15*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
19c4762a1bSJed Brown   for (i=0; i<m; i++) {
20c4762a1bSJed Brown     for (j=0; j<n; j++) {
21c4762a1bSJed Brown       v = 0.0 - 1.0*PETSC_i;
22*5f80ce2aSJacob Faibussowitsch       if (i>j && i-j<2)   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
23c4762a1bSJed Brown     }
24c4762a1bSJed Brown   }
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
27c4762a1bSJed Brown 
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHermitianTranspose(C, &C_htransposed));
29c4762a1bSJed Brown 
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_SELF));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C_htransposed,MAT_COPY_VALUES,&Cht));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(Cht,PETSC_VIEWER_STDOUT_SELF));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C_htransposed,MAT_DO_NOT_COPY_VALUES,&C_empty));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C_empty,PETSC_VIEWER_STDOUT_SELF));
35c4762a1bSJed Brown 
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C_htransposed));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Cht));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C_empty));
40c4762a1bSJed Brown   ierr = PetscFinalize();
41c4762a1bSJed Brown   return ierr;
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*TEST
45c4762a1bSJed Brown 
46c4762a1bSJed Brown    build:
47c4762a1bSJed Brown      requires: complex
48c4762a1bSJed Brown 
49c4762a1bSJed Brown    test:
50c4762a1bSJed Brown      output_file: output/ex175.out
51c4762a1bSJed Brown 
52c4762a1bSJed Brown TEST*/
53