xref: /petsc/src/mat/tests/ex253.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
16048efd1SBarry Smith static char help[] = "Tests MatMultHermitianTranspose() for real numbers.\n\n";
26048efd1SBarry Smith #include <petsc.h>
36048efd1SBarry Smith 
46048efd1SBarry Smith int main(int argc, char **args)
56048efd1SBarry Smith {
66048efd1SBarry Smith   Mat            A, AHT;
76048efd1SBarry Smith   Vec            x, y;
86048efd1SBarry Smith   PetscRandom    rand;
96048efd1SBarry Smith   PetscErrorCode ierr;
106048efd1SBarry Smith 
116048efd1SBarry Smith   ierr = PetscInitialize(&argc, &args, (char*)0, help); if (ierr) return ierr;
126048efd1SBarry Smith 
13*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
14*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A, MATAIJ));
15*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 10, 10));
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
186048efd1SBarry Smith 
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(A, 0, 0, 1.0, INSERT_VALUES));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
226048efd1SBarry Smith 
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHermitianTranspose(A, &AHT));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(AHT, &x, &y));
256048efd1SBarry Smith 
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(y, rand));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
306048efd1SBarry Smith 
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultHermitianTranspose(AHT, y, x));
326048efd1SBarry Smith 
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AHT));
376048efd1SBarry Smith   ierr = PetscFinalize();
386048efd1SBarry Smith   return ierr;
396048efd1SBarry Smith }
406048efd1SBarry Smith 
416048efd1SBarry Smith /*TEST
426048efd1SBarry Smith 
436048efd1SBarry Smith    test:
446048efd1SBarry Smith 
456048efd1SBarry Smith TEST*/
46