1c4762a1bSJed Brown static char help[] = "Test MatMultHermitianTranspose() and MatMultHermitianTransposeAdd().\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7e573050aSPierre Jolivet Mat A, B, C; 8c4762a1bSJed Brown Vec x, y, ys; 9c4762a1bSJed Brown PetscInt i, j; 10c4762a1bSJed Brown PetscScalar v; 11e573050aSPierre Jolivet PetscBool flg; 12c4762a1bSJed Brown 13327415f7SBarry Smith PetscFunctionBeginUser; 149566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 159566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 169566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 179566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 189566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 199566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 20c4762a1bSJed Brown 219371c9d4SSatish Balay i = 0; 229371c9d4SSatish Balay j = 0; 239371c9d4SSatish Balay v = 2.0; 249566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 259371c9d4SSatish Balay i = 0; 269371c9d4SSatish Balay j = 1; 279371c9d4SSatish Balay v = 3.0 + 4.0 * PETSC_i; 289566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 299371c9d4SSatish Balay i = 1; 309371c9d4SSatish Balay j = 0; 319371c9d4SSatish Balay v = 5.0 + 6.0 * PETSC_i; 329566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 339371c9d4SSatish Balay i = 1; 349371c9d4SSatish Balay j = 1; 359371c9d4SSatish Balay v = 7.0 + 8.0 * PETSC_i; 369566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Create vectors */ 419566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &y)); 429566063dSJacob Faibussowitsch PetscCall(VecSetSizes(y, PETSC_DECIDE, 2)); 439566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(y)); 449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ys)); 459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &x)); 46c4762a1bSJed Brown 479371c9d4SSatish Balay i = 0; 489371c9d4SSatish Balay v = 10.0 + 11.0 * PETSC_i; 499566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES)); 509371c9d4SSatish Balay i = 1; 519371c9d4SSatish Balay v = 100.0 + 120.0 * PETSC_i; 529566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES)); 539566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 549566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(A, x, y)); 579566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 589566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(A, x, y, ys)); 599566063dSJacob Faibussowitsch PetscCall(VecView(ys, PETSC_VIEWER_STDOUT_WORLD)); 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &B)); 629566063dSJacob Faibussowitsch PetscCall(MatCreateHermitianTranspose(A, &C)); 639566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeEqual(B, C, 4, &flg)); 6428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "B^Hx != C^Hx"); 659566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAddEqual(B, C, 4, &flg)); 6628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "y+B^Hx != y+C^Hx"); 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 69e573050aSPierre Jolivet 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 71c4762a1bSJed Brown 729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ys)); 759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 76b122ec5aSJacob Faibussowitsch return 0; 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /*TEST 80c4762a1bSJed Brown 81c4762a1bSJed Brown build: 82c4762a1bSJed Brown requires: complex 83*459e8d23SBlanca Mellado Pinto testset: 84c4762a1bSJed Brown test: 85*459e8d23SBlanca Mellado Pinto suffix: 1 86*459e8d23SBlanca Mellado Pinto args: -mat_type {{aij dense}} 87c4762a1bSJed Brown test: 88c4762a1bSJed Brown suffix: 2 89*459e8d23SBlanca Mellado Pinto args: -mat_type {{aij dense}} 90*459e8d23SBlanca Mellado Pinto diff_args: -j 91c4762a1bSJed Brown nsize: 2 92c4762a1bSJed Brown 93c4762a1bSJed Brown TEST*/ 94