Actual source code: ex19.c
petsc-3.4.5 2014-06-29
2: static char help[] = "Tests reusing MPI parallel matrices and MatGetValues().\n\
3: To test the parallel matrix assembly, this example intentionally lays out\n\
4: the matrix across processors differently from the way it is assembled.\n\
5: This example uses bilinear elements on the unit square. Input arguments are:\n\
6: -m <size> : problem size\n\n";
8: #include <petscmat.h>
12: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
13: {
15: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
16: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
17: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
18: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
19: return(0);
20: }
24: int main(int argc,char **args)
25: {
26: Mat C;
27: Vec u,b;
29: PetscMPIInt size,rank;
30: PetscInt i,m = 5,N,start,end,M,idx[4];
31: PetscInt j,nrsub,ncsub,*rsub,*csub,mystart,myend;
32: PetscBool flg;
33: PetscScalar one = 1.0,Ke[16],*vals;
34: PetscReal h,norm;
36: PetscInitialize(&argc,&args,(char*)0,help);
37: PetscOptionsGetInt(NULL,"-m",&m,NULL);
39: N = (m+1)*(m+1); /* dimension of matrix */
40: M = m*m; /* number of elements */
41: h = 1.0/m; /* mesh width */
42: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
43: MPI_Comm_size(PETSC_COMM_WORLD,&size);
45: /* Create stiffness matrix */
46: MatCreate(PETSC_COMM_WORLD,&C);
47: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
48: MatSetFromOptions(C);
49: MatSetUp(C);
51: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
52: end = start + M/size + ((M%size) > rank);
54: /* Form the element stiffness for the Laplacian */
55: FormElementStiffness(h*h,Ke);
56: for (i=start; i<end; i++) {
57: /* location of lower left corner of element */
58: /* node numbers for the four corners of element */
59: idx[0] = (m+1)*(i/m) + (i % m);
60: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
61: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
62: }
63: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
64: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
66: /* Assemble the matrix again */
67: MatZeroEntries(C);
69: for (i=start; i<end; i++) {
70: /* location of lower left corner of element */
71: /* node numbers for the four corners of element */
72: idx[0] = (m+1)*(i/m) + (i % m);
73: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
74: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
75: }
76: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
79: /* Create test vectors */
80: VecCreate(PETSC_COMM_WORLD,&u);
81: VecSetSizes(u,PETSC_DECIDE,N);
82: VecSetFromOptions(u);
83: VecDuplicate(u,&b);
84: VecSet(u,one);
86: /* Check error */
87: MatMult(C,u,b);
88: VecNorm(b,NORM_2,&norm);
89: if (norm > 1.e-10 || norm < -1.e-10) {
90: PetscPrintf(PETSC_COMM_WORLD,"Norm of error b %G should be near 0\n",norm);
91: }
93: /* Now test MatGetValues() */
94: PetscOptionsHasName(NULL,"-get_values",&flg);
95: if (flg) {
96: MatGetOwnershipRange(C,&mystart,&myend);
97: nrsub = myend - mystart; ncsub = 4;
98: PetscMalloc(nrsub*ncsub*sizeof(PetscScalar),&vals);
99: PetscMalloc(nrsub*sizeof(PetscInt),&rsub);
100: PetscMalloc(ncsub*sizeof(PetscInt),&csub);
101: for (i=myend-1; i>=mystart; i--) rsub[myend-i-1] = i;
102: for (i=0; i<ncsub; i++) csub[i] = 2*(ncsub-i) + mystart;
103: MatGetValues(C,nrsub,rsub,ncsub,csub,vals);
104: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
105: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"processor number %d: start=%D, end=%D, mystart=%D, myend=%D\n",rank,start,end,mystart,myend);
106: for (i=0; i<nrsub; i++) {
107: for (j=0; j<ncsub; j++) {
108: if (PetscImaginaryPart(vals[i*ncsub+j]) != 0.0) {
109: PetscSynchronizedPrintf(PETSC_COMM_WORLD," C[%D, %D] = %G + %G i\n",rsub[i],csub[j],PetscRealPart(vals[i*ncsub+j]),PetscImaginaryPart(vals[i*ncsub+j]));
110: } else {
111: PetscSynchronizedPrintf(PETSC_COMM_WORLD," C[%D, %D] = %G\n",rsub[i],csub[j],PetscRealPart(vals[i*ncsub+j]));
112: }
113: }
114: }
115: PetscSynchronizedFlush(PETSC_COMM_WORLD);
116: PetscFree(rsub);
117: PetscFree(csub);
118: PetscFree(vals);
119: }
121: /* Free data structures */
122: VecDestroy(&u);
123: VecDestroy(&b);
124: MatDestroy(&C);
125: PetscFinalize();
126: return 0;
127: }