xref: /petsc/src/dm/tests/ex30.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests DMSLICED operations\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmsliced.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown int main(int argc,char *argv[])
6*c4762a1bSJed Brown {
7*c4762a1bSJed Brown   char           mat_type[256] = MATAIJ; /* default matrix type */
8*c4762a1bSJed Brown   PetscErrorCode ierr;
9*c4762a1bSJed Brown   MPI_Comm       comm;
10*c4762a1bSJed Brown   PetscMPIInt    rank,size;
11*c4762a1bSJed Brown   DM             slice;
12*c4762a1bSJed Brown   PetscInt       i,bs=1,N=5,n,m,rstart,ghosts[2],*d_nnz,*o_nnz,dfill[4]={1,0,0,1},ofill[4]={1,1,1,1};
13*c4762a1bSJed Brown   PetscReal      alpha   =1,K=1,rho0=1,u0=0,sigma=0.2;
14*c4762a1bSJed Brown   PetscBool      useblock=PETSC_TRUE;
15*c4762a1bSJed Brown   PetscScalar    *xx;
16*c4762a1bSJed Brown   Mat            A;
17*c4762a1bSJed Brown   Vec            x,b,lf;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
20*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
21*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,0,"Options for DMSliced test",0);CHKERRQ(ierr);
25*c4762a1bSJed Brown   {
26*c4762a1bSJed Brown     ierr = PetscOptionsInt("-n","Global number of nodes","",N,&N,NULL);CHKERRQ(ierr);
27*c4762a1bSJed Brown     ierr = PetscOptionsInt("-bs","Block size (1 or 2)","",bs,&bs,NULL);CHKERRQ(ierr);
28*c4762a1bSJed Brown     if (bs != 1) {
29*c4762a1bSJed Brown       if (bs != 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Block size must be 1 or 2");
30*c4762a1bSJed Brown       ierr = PetscOptionsReal("-alpha","Inverse time step for wave operator","",alpha,&alpha,NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown       ierr = PetscOptionsReal("-K","Bulk modulus of compressibility","",K,&K,NULL);CHKERRQ(ierr);
32*c4762a1bSJed Brown       ierr = PetscOptionsReal("-rho0","Reference density","",rho0,&rho0,NULL);CHKERRQ(ierr);
33*c4762a1bSJed Brown       ierr = PetscOptionsReal("-u0","Reference velocity","",u0,&u0,NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown       ierr = PetscOptionsReal("-sigma","Width of Gaussian density perturbation","",sigma,&sigma,NULL);CHKERRQ(ierr);
35*c4762a1bSJed Brown       ierr = PetscOptionsBool("-block","Use block matrix assembly","",useblock,&useblock,NULL);CHKERRQ(ierr);
36*c4762a1bSJed Brown     }
37*c4762a1bSJed Brown     ierr = PetscOptionsString("-sliced_mat_type","Matrix type to use (aij or baij)","",mat_type,mat_type,sizeof(mat_type),NULL);CHKERRQ(ierr);
38*c4762a1bSJed Brown   }
39*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown   /* Split ownership, set up periodic grid in 1D */
42*c4762a1bSJed Brown   n         = PETSC_DECIDE;
43*c4762a1bSJed Brown   ierr      = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr);
44*c4762a1bSJed Brown   rstart    = 0;
45*c4762a1bSJed Brown   ierr      = MPI_Scan(&n,&rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
46*c4762a1bSJed Brown   rstart   -= n;
47*c4762a1bSJed Brown   ghosts[0] = (N+rstart-1)%N;
48*c4762a1bSJed Brown   ghosts[1] = (rstart+n)%N;
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   ierr = PetscMalloc2(n,&d_nnz,n,&o_nnz);CHKERRQ(ierr);
51*c4762a1bSJed Brown   for (i=0; i<n; i++) {
52*c4762a1bSJed Brown     if (size > 1 && (i==0 || i==n-1)) {
53*c4762a1bSJed Brown       d_nnz[i] = 2;
54*c4762a1bSJed Brown       o_nnz[i] = 1;
55*c4762a1bSJed Brown     } else {
56*c4762a1bSJed Brown       d_nnz[i] = 3;
57*c4762a1bSJed Brown       o_nnz[i] = 0;
58*c4762a1bSJed Brown     }
59*c4762a1bSJed Brown   }
60*c4762a1bSJed Brown   ierr = DMSlicedCreate(comm,bs,n,2,ghosts,d_nnz,o_nnz,&slice);CHKERRQ(ierr); /* Currently does not copy X_nnz so we can't free them until after DMSlicedGetMatrix */
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   if (!useblock) {ierr = DMSlicedSetBlockFills(slice,dfill,ofill);CHKERRQ(ierr);} /* Irrelevant for baij formats */
63*c4762a1bSJed Brown   ierr = DMSetMatType(slice,mat_type);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = DMCreateMatrix(slice,&A);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(slice,&x);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = VecDuplicate(x,&b);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   ierr = VecGhostGetLocalForm(x,&lf);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = VecGetSize(lf,&m);CHKERRQ(ierr);
73*c4762a1bSJed Brown   if (m != (n+2)*bs) SETERRQ2(PETSC_COMM_SELF,1,"size of local form %D, expected %D",m,(n+2)*bs);
74*c4762a1bSJed Brown   ierr = VecGetArray(lf,&xx);CHKERRQ(ierr);
75*c4762a1bSJed Brown   for (i=0; i<n; i++) {
76*c4762a1bSJed Brown     PetscInt        row[2],col[9],im,ip;
77*c4762a1bSJed Brown     PetscScalar     v[12];
78*c4762a1bSJed Brown     const PetscReal xref = 2.0*(rstart+i)/N - 1; /* [-1,1] */
79*c4762a1bSJed Brown     const PetscReal h    = 1.0/N;                /* grid spacing */
80*c4762a1bSJed Brown     im = (i==0) ? n : i-1;
81*c4762a1bSJed Brown     ip = (i==n-1) ? n+1 : i+1;
82*c4762a1bSJed Brown     switch (bs) {
83*c4762a1bSJed Brown     case 1:                     /* Laplacian with periodic boundaries */
84*c4762a1bSJed Brown       col[0] = im;         col[1] = i;        col[2] = ip;
85*c4762a1bSJed Brown       v[0]   = -h;           v[1] = 2*h;        v[2] = -h;
86*c4762a1bSJed Brown       ierr   = MatSetValuesLocal(A,1,&i,3,col,v,INSERT_VALUES);CHKERRQ(ierr);
87*c4762a1bSJed Brown       xx[i]  = PetscSinReal(xref*PETSC_PI);
88*c4762a1bSJed Brown       break;
89*c4762a1bSJed Brown     case 2:                     /* Linear acoustic wave operator in variables [rho, u], central differences, periodic, timestep 1/alpha */
90*c4762a1bSJed Brown       v[0] = -0.5*u0;   v[1] = -0.5*K;      v[2] = alpha; v[3] = 0;       v[4] = 0.5*u0;    v[5] = 0.5*K;
91*c4762a1bSJed Brown       v[6] = -0.5/rho0; v[7] = -0.5*u0;     v[8] = 0;     v[9] = alpha;   v[10] = 0.5/rho0; v[11] = 0.5*u0;
92*c4762a1bSJed Brown       if (useblock) {
93*c4762a1bSJed Brown         row[0] = i; col[0] = im; col[1] = i; col[2] = ip;
94*c4762a1bSJed Brown         ierr   = MatSetValuesBlockedLocal(A,1,row,3,col,v,INSERT_VALUES);CHKERRQ(ierr);
95*c4762a1bSJed Brown       } else {
96*c4762a1bSJed Brown         row[0] = 2*i; row[1] = 2*i+1;
97*c4762a1bSJed Brown         col[0] = 2*im; col[1] = 2*im+1; col[2] = 2*i; col[3] = 2*ip; col[4] = 2*ip+1;
98*c4762a1bSJed Brown         v[3]   = v[4]; v[4] = v[5];                                                     /* pack values in first row */
99*c4762a1bSJed Brown         ierr   = MatSetValuesLocal(A,1,row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
100*c4762a1bSJed Brown         col[2] = 2*i+1;
101*c4762a1bSJed Brown         v[8]   = v[9]; v[9] = v[10]; v[10] = v[11];                                     /* pack values in second row */
102*c4762a1bSJed Brown         ierr   = MatSetValuesLocal(A,1,row+1,5,col,v+6,INSERT_VALUES);CHKERRQ(ierr);
103*c4762a1bSJed Brown       }
104*c4762a1bSJed Brown       /* Set current state (gaussian density perturbation) */
105*c4762a1bSJed Brown       xx[2*i]   = 0.2*PetscExpReal(-PetscSqr(xref)/(2*PetscSqr(sigma)));
106*c4762a1bSJed Brown       xx[2*i+1] = 0;
107*c4762a1bSJed Brown       break;
108*c4762a1bSJed Brown     default: SETERRQ1(PETSC_COMM_SELF,1,"not implemented for block size %D",bs);
109*c4762a1bSJed Brown     }
110*c4762a1bSJed Brown   }
111*c4762a1bSJed Brown   ierr = VecRestoreArray(lf,&xx);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = VecGhostRestoreLocalForm(x,&lf);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   ierr = MatMult(A,x,b);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = VecView(b,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown   /* Update the ghosted values, view the result on rank 0. */
122*c4762a1bSJed Brown   ierr = VecGhostUpdateBegin(b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = VecGhostUpdateEnd(b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124*c4762a1bSJed Brown   if (!rank) {
125*c4762a1bSJed Brown     ierr = VecGhostGetLocalForm(b,&lf);CHKERRQ(ierr);
126*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"Local form of b on rank 0, last two nodes are ghost nodes\n");CHKERRQ(ierr);
127*c4762a1bSJed Brown     ierr = VecView(lf,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
128*c4762a1bSJed Brown     ierr = VecGhostRestoreLocalForm(b,&lf);CHKERRQ(ierr);
129*c4762a1bSJed Brown   }
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown   ierr = DMDestroy(&slice);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = PetscFinalize();
136*c4762a1bSJed Brown   return ierr;
137*c4762a1bSJed Brown }
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown /*TEST
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown    test:
143*c4762a1bSJed Brown       nsize: 2
144*c4762a1bSJed Brown       args: -bs 2 -block 0 -sliced_mat_type baij -alpha 10 -u0 0.1
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown    test:
147*c4762a1bSJed Brown       suffix: 2
148*c4762a1bSJed Brown       nsize: 2
149*c4762a1bSJed Brown       args: -bs 2 -block 1 -sliced_mat_type aij -alpha 10 -u0 0.1
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown    test:
152*c4762a1bSJed Brown       suffix: 3
153*c4762a1bSJed Brown       nsize: 2
154*c4762a1bSJed Brown       args: -bs 2 -block 0 -sliced_mat_type aij -alpha 10 -u0 0.1
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown TEST*/
157