Actual source code: ex175.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: static char help[] = "Tests MatCreateHermitianTranspose().\n\n";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  8:   Mat            C,C_htransposed,Cht,C_empty;
  9:   PetscInt       i,j,m = 10,n = 10;
 11:   PetscScalar    v;

 13:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 14:   /* Create a complex non-hermitian matrix */
 15:   MatCreate(PETSC_COMM_SELF,&C);
 16:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n);
 17:   MatSetFromOptions(C);
 18:   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)   {MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);}
 23:     }
 24:   }
 25:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 26:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 28:   MatCreateHermitianTranspose(C, &C_htransposed);

 30:   MatView(C,PETSC_VIEWER_STDOUT_SELF);
 31:   MatDuplicate(C_htransposed,MAT_COPY_VALUES,&Cht);
 32:   MatView(Cht,PETSC_VIEWER_STDOUT_SELF);
 33:   MatDuplicate(C_htransposed,MAT_DO_NOT_COPY_VALUES,&C_empty);
 34:   MatView(C_empty,PETSC_VIEWER_STDOUT_SELF);

 36:   MatDestroy(&C);
 37:   MatDestroy(&C_htransposed);
 38:   MatDestroy(&Cht);
 39:   MatDestroy(&C_empty);
 40:   PetscFinalize();
 41:   return ierr;
 42: }



 46: /*TEST

 48:    build:
 49:      requires: complex

 51:    test:
 52:      output_file: output/ex175.out

 54: TEST*/