xref: /petsc/src/dm/tests/ex22.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 static char help[] = "Tests MatSetValuesBlockedStencil() in 3d.\n\n";
3 
4 #include <petscmat.h>
5 #include <petscdm.h>
6 #include <petscdmda.h>
7 
8 int main(int argc,char **argv)
9 {
10   PetscInt        M = 3,N = 4,P = 2,s = 1,w = 2,i, m = PETSC_DECIDE,n = PETSC_DECIDE,p = PETSC_DECIDE;
11   DM              da;
12   Mat             mat;
13   DMDAStencilType stencil_type = DMDA_STENCIL_BOX;
14   PetscBool       flg          = PETSC_FALSE;
15   MatStencil      idx[2],idy[2];
16   PetscScalar     *values;
17 
18   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
19   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
20   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
21   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-P",&P,NULL));
22   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
23   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
24   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
25   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-s",&s,NULL));
26   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-w",&w,NULL));
27   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-star",&flg,NULL));
28   if (flg) stencil_type =  DMDA_STENCIL_STAR;
29 
30   /* Create distributed array and get vectors */
31   CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,stencil_type,M,N,P,m,n,p,w,s,0,0,0,&da));
32   CHKERRQ(DMSetFromOptions(da));
33   CHKERRQ(DMSetUp(da));
34   CHKERRQ(DMSetMatType(da,MATMPIBAIJ));
35   CHKERRQ(DMCreateMatrix(da,&mat));
36 
37   idx[0].i = 1;   idx[0].j = 1; idx[0].k = 0;
38   idx[1].i = 2;   idx[1].j = 1; idx[1].k = 0;
39   idy[0].i = 1;   idy[0].j = 2; idy[0].k = 0;
40   idy[1].i = 2;   idy[1].j = 2; idy[1].k = 0;
41   CHKERRQ(PetscMalloc1(2*2*w*w,&values));
42   for (i=0; i<2*2*w*w; i++) values[i] = i;
43   CHKERRQ(MatSetValuesBlockedStencil(mat,2,idx,2,idy,values,INSERT_VALUES));
44   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
45   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
46 
47   /* Free memory */
48   CHKERRQ(PetscFree(values));
49   CHKERRQ(MatDestroy(&mat));
50   CHKERRQ(DMDestroy(&da));
51   CHKERRQ(PetscFinalize());
52   return 0;
53 }
54