xref: /petsc/src/mat/tests/ex253.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
10*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &args, (char*)0, help));
116048efd1SBarry Smith 
125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A, MATAIJ));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 10, 10));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
176048efd1SBarry Smith 
185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(A, 0, 0, 1.0, INSERT_VALUES));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
216048efd1SBarry Smith 
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHermitianTranspose(A, &AHT));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(AHT, &x, &y));
246048efd1SBarry Smith 
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(y, rand));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
296048efd1SBarry Smith 
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultHermitianTranspose(AHT, y, x));
316048efd1SBarry Smith 
325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AHT));
36*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
37*b122ec5aSJacob Faibussowitsch   return 0;
386048efd1SBarry Smith }
396048efd1SBarry Smith 
406048efd1SBarry Smith /*TEST
416048efd1SBarry Smith 
426048efd1SBarry Smith    test:
436048efd1SBarry Smith 
446048efd1SBarry Smith TEST*/
45