1c4762a1bSJed Brown static char help[] = "Tests MatCreateHermitianTranspose().\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7a1f56445SPierre Jolivet Mat C, C_htransposed, Cht, C_empty, Cht_empty; 8c4762a1bSJed Brown PetscInt i, j, m = 10, n = 10; 9c4762a1bSJed Brown PetscScalar v; 10c4762a1bSJed Brown 11327415f7SBarry Smith PetscFunctionBeginUser; 12*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 13c4762a1bSJed Brown /* Create a complex non-hermitian matrix */ 149566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &C)); 159566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m, n)); 169566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 179566063dSJacob Faibussowitsch PetscCall(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; 219566063dSJacob Faibussowitsch if (i > j && i - j < 2) PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES)); 22c4762a1bSJed Brown } 23c4762a1bSJed Brown } 249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 26c4762a1bSJed Brown 279566063dSJacob Faibussowitsch PetscCall(MatCreateHermitianTranspose(C, &C_htransposed)); 28c4762a1bSJed Brown 299566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C_htransposed, MAT_COPY_VALUES, &Cht)); 31a1f56445SPierre Jolivet PetscCall(MatConvert(Cht, MATAIJ, MAT_INPLACE_MATRIX, &Cht)); 329566063dSJacob Faibussowitsch PetscCall(MatView(Cht, PETSC_VIEWER_STDOUT_SELF)); 339566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C_htransposed, MAT_DO_NOT_COPY_VALUES, &C_empty)); 34a1f56445SPierre Jolivet PetscCall(MatHermitianTransposeGetMat(C_empty, &Cht_empty)); 35a1f56445SPierre Jolivet PetscCall(MatHermitianTranspose(Cht_empty, MAT_INPLACE_MATRIX, &Cht_empty)); 36a1f56445SPierre Jolivet PetscCall(MatView(Cht_empty, PETSC_VIEWER_STDOUT_SELF)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_htransposed)); 409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cht)); 419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_empty)); 429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 43b122ec5aSJacob Faibussowitsch return 0; 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown /*TEST 47c4762a1bSJed Brown 48c4762a1bSJed Brown build: 49c4762a1bSJed Brown requires: complex 50c4762a1bSJed Brown 51c4762a1bSJed Brown test: 52c4762a1bSJed Brown output_file: output/ex175.out 53c4762a1bSJed Brown 54c4762a1bSJed Brown TEST*/ 55