xref: /petsc/src/dm/tests/ex5.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown /*
2*c4762a1bSJed Brown   Test DMCreateMatrix() for structure_only
3*c4762a1bSJed Brown */
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #include <petscdmda.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown int main(int argc, char *argv[])
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscErrorCode ierr;
10*c4762a1bSJed Brown   PetscInt       nx=6,ny=6,nz=6,dim=1,dof=2;
11*c4762a1bSJed Brown   DM             da;
12*c4762a1bSJed Brown   Mat            A;
13*c4762a1bSJed Brown   PetscBool      struct_only=PETSC_TRUE;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,NULL);if (ierr) return ierr;
16*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr);
17*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);CHKERRQ(ierr);
18*c4762a1bSJed Brown   switch (dim) {
19*c4762a1bSJed Brown   case 1:
20*c4762a1bSJed Brown     ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,nx,dof,1,NULL,&da);CHKERRQ(ierr);
21*c4762a1bSJed Brown     break;
22*c4762a1bSJed Brown   case 2:
23*c4762a1bSJed Brown     ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,nx,ny,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&da);CHKERRQ(ierr);
24*c4762a1bSJed Brown     break;
25*c4762a1bSJed Brown   default:
26*c4762a1bSJed Brown     ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nx,ny,nz,
27*c4762a1bSJed Brown                       PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,2,NULL,NULL,NULL,&da);CHKERRQ(ierr);
28*c4762a1bSJed Brown   }
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = DMView(da,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-struct_only",&struct_only,NULL);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = DMSetMatrixStructureOnly(da,struct_only);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr);
37*c4762a1bSJed Brown   /* Set da->structure_only to default PETSC_FALSE in case da is being used to create new matrices */
38*c4762a1bSJed Brown   ierr = DMSetMatrixStructureOnly(da,PETSC_FALSE);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscFinalize();
44*c4762a1bSJed Brown   return ierr;
45*c4762a1bSJed Brown }
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown /*TEST
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown    test:
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown    test:
53*c4762a1bSJed Brown       suffix: 2
54*c4762a1bSJed Brown       args: -dm_mat_type baij -dim 2
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown TEST*/
57