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; 14*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, 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 */ 4195571869SBlanca Mellado Pinto PetscCall(MatCreateVecs(A, &x, &y)); 429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ys)); 43c4762a1bSJed Brown 449371c9d4SSatish Balay i = 0; 459371c9d4SSatish Balay v = 10.0 + 11.0 * PETSC_i; 469566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES)); 479371c9d4SSatish Balay i = 1; 489371c9d4SSatish Balay v = 100.0 + 120.0 * PETSC_i; 499566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES)); 509566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 519566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(A, x, y)); 549566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 559566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(A, x, y, ys)); 569566063dSJacob Faibussowitsch PetscCall(VecView(ys, PETSC_VIEWER_STDOUT_WORLD)); 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &B)); 599566063dSJacob Faibussowitsch PetscCall(MatCreateHermitianTranspose(A, &C)); 609566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeEqual(B, C, 4, &flg)); 6128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "B^Hx != C^Hx"); 629566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAddEqual(B, C, 4, &flg)); 6328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "y+B^Hx != y+C^Hx"); 649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 66e573050aSPierre Jolivet 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ys)); 729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 73b122ec5aSJacob Faibussowitsch return 0; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown /*TEST 77c4762a1bSJed Brown 78c4762a1bSJed Brown build: 79c4762a1bSJed Brown requires: complex 805e4d33a3SBlanca Mellado Pinto 81459e8d23SBlanca Mellado Pinto testset: 825e4d33a3SBlanca Mellado Pinto output_file: output/ex197_1.out 83c4762a1bSJed Brown test: 84459e8d23SBlanca Mellado Pinto suffix: 1 85459e8d23SBlanca Mellado Pinto args: -mat_type {{aij dense}} 86c4762a1bSJed Brown test: 8795571869SBlanca Mellado Pinto suffix: 1_cuda 8895571869SBlanca Mellado Pinto requires: cuda 8995571869SBlanca Mellado Pinto args: -mat_type densecuda 9095571869SBlanca Mellado Pinto filter: sed -e 's/seqcuda/seq/' 915e4d33a3SBlanca Mellado Pinto 925e4d33a3SBlanca Mellado Pinto testset: 935e4d33a3SBlanca Mellado Pinto output_file: output/ex197_2.out 945e4d33a3SBlanca Mellado Pinto nsize: 2 9595571869SBlanca Mellado Pinto test: 96c4762a1bSJed Brown suffix: 2 97459e8d23SBlanca Mellado Pinto args: -mat_type {{aij dense}} 985e4d33a3SBlanca Mellado Pinto test: 995e4d33a3SBlanca Mellado Pinto suffix: 2_scalapack 1005e4d33a3SBlanca Mellado Pinto requires: scalapack 1015e4d33a3SBlanca Mellado Pinto args: -mat_type scalapack 102c4762a1bSJed Brown 103c4762a1bSJed Brown TEST*/ 104