xref: /petsc/src/dm/tests/ex31.c (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1 static char help[] = "Tests MAIJ matrix for large DOF\n\n";
2 
3 #include <petscdm.h>
4 #include <petscdmda.h>
5 
6 int main(int argc, char *argv[])
7 {
8   Mat M;
9   Vec x, y;
10   DM  da, daf;
11 
12   PetscFunctionBeginUser;
13   PetscCall(PetscInitialize(&argc, &argv, 0, help));
14   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 5, PETSC_DECIDE, PETSC_DECIDE, 41, 1, 0, 0, &da));
15   PetscCall(DMSetFromOptions(da));
16   PetscCall(DMSetUp(da));
17   PetscCall(DMRefine(da, PETSC_COMM_WORLD, &daf));
18   PetscCall(DMCreateInterpolation(da, daf, &M, NULL));
19   PetscCall(DMCreateGlobalVector(da, &x));
20   PetscCall(DMCreateGlobalVector(daf, &y));
21 
22   PetscCall(MatMult(M, x, y));
23   PetscCall(MatMultTranspose(M, y, x));
24   PetscCall(DMDestroy(&da));
25   PetscCall(DMDestroy(&daf));
26   PetscCall(VecDestroy(&x));
27   PetscCall(VecDestroy(&y));
28   PetscCall(MatDestroy(&M));
29   PetscCall(PetscFinalize());
30   return 0;
31 }
32 
33 /*TEST
34 
35    test:
36       nsize: 2
37       output_file: output/empty.out
38 
39 TEST*/
40