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