xref: /petsc/src/mat/tests/ex253.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
16048efd1SBarry Smith static char help[] = "Tests MatMultHermitianTranspose() for real numbers.\n\n";
26048efd1SBarry Smith #include <petsc.h>
36048efd1SBarry Smith 
main(int argc,char ** args)4d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
5d71ae5a4SJacob Faibussowitsch {
66048efd1SBarry Smith   Mat         A, AHT;
76048efd1SBarry Smith   Vec         x, y;
86048efd1SBarry Smith   PetscRandom rand;
96048efd1SBarry Smith 
10327415f7SBarry Smith   PetscFunctionBeginUser;
11c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
126048efd1SBarry Smith 
139566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
149566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
159566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 10, 10));
169566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
179566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
186048efd1SBarry Smith 
199566063dSJacob Faibussowitsch   PetscCall(MatSetValue(A, 0, 0, 1.0, INSERT_VALUES));
209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
226048efd1SBarry Smith 
239566063dSJacob Faibussowitsch   PetscCall(MatCreateHermitianTranspose(A, &AHT));
249566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(AHT, &x, &y));
256048efd1SBarry Smith 
269566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
289566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(y, rand));
299566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
306048efd1SBarry Smith 
319566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(AHT, y, x));
326048efd1SBarry Smith 
339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AHT));
379566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
38b122ec5aSJacob Faibussowitsch   return 0;
396048efd1SBarry Smith }
406048efd1SBarry Smith 
416048efd1SBarry Smith /*TEST
426048efd1SBarry Smith 
436048efd1SBarry Smith    test:
44*3886731fSPierre Jolivet      output_file: output/empty.out
456048efd1SBarry Smith 
466048efd1SBarry Smith TEST*/
47