xref: /petsc/src/mat/tests/ex247.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 
2 static char help[] = "Tests MATCENTERING matrix type.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **argv)
7 {
8   PetscErrorCode ierr;
9   PetscInt       n;
10   Mat            C;
11   Vec            x,y;
12   PetscReal      norm;
13   PetscMPIInt    size;
14 
15   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
17 
18   /* Create a parallel vector with 10*size total entries, and fill it with 1s. */
19   n = 10*size;
20   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
21   ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr);
22   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
23   ierr = VecSet(x,1.0);CHKERRQ(ierr);
24 
25   /* Create a corresponding n x n centering matrix and use it to create a mean-centered y = C * x. */
26   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
27   ierr = MatCreateCentering(PETSC_COMM_WORLD,PETSC_DECIDE,n,&C);CHKERRQ(ierr);
28   ierr = MatMult(C,x,y);CHKERRQ(ierr);
29 
30   /* Verify that the centered vector y has norm 0. */
31   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
32   ierr = PetscPrintf(PETSC_COMM_WORLD,"Vector norm after MatMult() with centering matrix applied to vector of ones is %f.\n",(double)norm);CHKERRQ(ierr);
33 
34   /* Now repeat, but using MatMultTranspose(). */
35   ierr = MatMultTranspose(C,x,y);CHKERRQ(ierr);
36   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
37   ierr = PetscPrintf(PETSC_COMM_WORLD,"Vector norm after MatMultTranspose() with centering matrix applied to vector of ones is %f.\n",(double)norm);CHKERRQ(ierr);
38 
39   /* Clean up. */
40   ierr = VecDestroy(&x);CHKERRQ(ierr);
41   ierr = VecDestroy(&y);CHKERRQ(ierr);
42   ierr = MatDestroy(&C);CHKERRQ(ierr);
43   ierr = PetscFinalize();
44   return ierr;
45 }
46 
47 /*TEST
48 
49     test:
50       suffix: 1
51       nsize: 1
52       output_file: output/ex247.out
53 
54     test:
55       suffix: 2
56       nsize: 2
57       output_file: output/ex247.out
58 
59 TEST*/
60