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