Actual source code: ex3.c


  2: static char help[] = "Bilinear elements on the unit square for Laplacian.  To test the parallel\n\
  3: matrix assembly, the matrix is intentionally laid out across processors\n\
  4: differently from the way it is assembled.  Input arguments are:\n\
  5:   -m <size> : problem size\n\n";

  7: /* Addendum: piggy-backing on this example to test KSPChebyshev methods */

  9: #include <petscksp.h>

 11: int FormElementStiffness(PetscReal H,PetscScalar *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: }
 20: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
 21: {
 23:   r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
 24:   return 0;
 25: }

 27: int main(int argc,char **args)
 28: {
 29:   Mat            C;
 30:   PetscMPIInt    rank,size;
 31:   PetscInt       i,m = 5,N,start,end,M,its;
 32:   PetscScalar    val,Ke[16],r[4];
 33:   PetscReal      x,y,h,norm;
 34:   PetscInt       idx[4],count,*rows;
 35:   Vec            u,ustar,b;
 36:   KSP            ksp;
 37:   PetscBool      viewkspest = PETSC_FALSE;

 39:   PetscInitialize(&argc,&args,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 41:   PetscOptionsGetBool(NULL,NULL,"-ksp_est_view",&viewkspest,NULL);
 42:   N    = (m+1)*(m+1); /* dimension of matrix */
 43:   M    = m*m; /* number of elements */
 44:   h    = 1.0/m;    /* mesh width */
 45:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 46:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 48:   /* Create stiffness matrix */
 49:   MatCreate(PETSC_COMM_WORLD,&C);
 50:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 51:   MatSetFromOptions(C);
 52:   MatSetUp(C);
 53:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 54:   end   = start + M/size + ((M%size) > rank);

 56:   /* Assemble matrix */
 57:   FormElementStiffness(h*h,Ke);   /* element stiffness for Laplacian */
 58:   for (i=start; i<end; i++) {
 59:     /* node numbers for the four corners of element */
 60:     idx[0] = (m+1)*(i/m) + (i % m);
 61:     idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 62:     MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
 63:   }
 64:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 65:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 67:   /* Create right-hand-side and solution vectors */
 68:   VecCreate(PETSC_COMM_WORLD,&u);
 69:   VecSetSizes(u,PETSC_DECIDE,N);
 70:   VecSetFromOptions(u);
 71:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 72:   VecDuplicate(u,&b);
 73:   PetscObjectSetName((PetscObject)b,"Right hand side");
 74:   VecDuplicate(b,&ustar);
 75:   VecSet(u,0.0);
 76:   VecSet(b,0.0);

 78:   /* Assemble right-hand-side vector */
 79:   for (i=start; i<end; i++) {
 80:     /* location of lower left corner of element */
 81:     x = h*(i % m); y = h*(i/m);
 82:     /* node numbers for the four corners of element */
 83:     idx[0] = (m+1)*(i/m) + (i % m);
 84:     idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 85:     FormElementRhs(x,y,h*h,r);
 86:     VecSetValues(b,4,idx,r,ADD_VALUES);
 87:   }
 88:   VecAssemblyBegin(b);
 89:   VecAssemblyEnd(b);

 91:   /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
 92:   PetscMalloc1(4*m,&rows);
 93:   for (i=0; i<m+1; i++) {
 94:     rows[i]          = i; /* bottom */
 95:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
 96:   }
 97:   count = m+1; /* left side */
 98:   for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;

100:   count = 2*m; /* left side */
101:   for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
102:   for (i=0; i<4*m; i++) {
103:     val  = h*(rows[i]/(m+1));
104:     VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
105:     VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
106:   }
107:   MatZeroRows(C,4*m,rows,1.0,0,0);

109:   PetscFree(rows);
110:   VecAssemblyBegin(u);
111:   VecAssemblyEnd(u);
112:   VecAssemblyBegin(b);
113:   VecAssemblyEnd(b);

115:   { Mat A;
116:     MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
117:     MatDestroy(&C);
118:     MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);
119:     MatDestroy(&A);
120:   }

122:   /* Solve linear system */
123:   KSPCreate(PETSC_COMM_WORLD,&ksp);
124:   KSPSetOperators(ksp,C,C);
125:   KSPSetFromOptions(ksp);
126:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
127:   KSPSolve(ksp,b,u);

129:   if (viewkspest) {
130:     KSP kspest;

132:     KSPChebyshevEstEigGetKSP(ksp,&kspest);
133:     if (kspest) KSPView(kspest,PETSC_VIEWER_STDOUT_WORLD);
134:   }

136:   /* Check error */
137:   VecGetOwnershipRange(ustar,&start,&end);
138:   for (i=start; i<end; i++) {
139:     val  = h*(i/(m+1));
140:     VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
141:   }
142:   VecAssemblyBegin(ustar);
143:   VecAssemblyEnd(ustar);
144:   VecAXPY(u,-1.0,ustar);
145:   VecNorm(u,NORM_2,&norm);
146:   KSPGetIterationNumber(ksp,&its);
147:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);

149:   /* Free work space */
150:   KSPDestroy(&ksp);
151:   VecDestroy(&ustar);
152:   VecDestroy(&u);
153:   VecDestroy(&b);
154:   MatDestroy(&C);
155:   PetscFinalize();
156:   return 0;
157: }

159: /*TEST

161:     test:
162:       args: -pc_type jacobi -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always

164:     test:
165:       suffix: 2
166:       nsize: 2
167:       args: -pc_type jacobi -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always

169:     test:
170:       suffix: 2_kokkos
171:       nsize: 2
172:       args: -pc_type jacobi -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always -mat_type aijkokkos -vec_type kokkos
173:       output_file: output/ex3_2.out
174:       requires: kokkos_kernels

176:     test:
177:       suffix: nocheby
178:       args: -ksp_est_view

180:     test:
181:       suffix: chebynoest
182:       args: -ksp_est_view -ksp_type chebyshev -ksp_chebyshev_eigenvalues 0.1,1.0

184:     test:
185:       suffix: chebyest
186:       args: -ksp_est_view -ksp_type chebyshev -ksp_chebyshev_esteig
187:       filter:  sed -e "s/Iterations 19/Iterations 20/g"

189: TEST*/