1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the multigrid code. The input parameters are:\n\ 3c4762a1bSJed Brown -x N Use a mesh in the x direction of N. \n\ 4c4762a1bSJed Brown -c N Use N V-cycles. \n\ 5c4762a1bSJed Brown -l N Use N Levels. \n\ 6c4762a1bSJed Brown -smooths N Use N pre smooths and N post smooths. \n\ 7c4762a1bSJed Brown -j Use Jacobi smoother. \n\ 8c4762a1bSJed Brown -a use additive multigrid \n\ 9c4762a1bSJed Brown -f use full multigrid (preconditioner variant) \n\ 10c4762a1bSJed Brown This example also demonstrates matrix-free methods\n\n"; 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown This is not a good example to understand the use of multigrid with PETSc. 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16c4762a1bSJed Brown #include <petscksp.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscErrorCode residual(Mat,Vec,Vec,Vec); 19c4762a1bSJed Brown PetscErrorCode gauss_seidel(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*); 20c4762a1bSJed Brown PetscErrorCode jacobi_smoother(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*); 21c4762a1bSJed Brown PetscErrorCode interpolate(Mat,Vec,Vec,Vec); 22c4762a1bSJed Brown PetscErrorCode restrct(Mat,Vec,Vec); 23c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt,Mat*); 24c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec); 25c4762a1bSJed Brown PetscErrorCode CalculateError(Vec,Vec,Vec,PetscReal*); 26c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt,Vec*); 27c4762a1bSJed Brown PetscErrorCode amult(Mat,Vec,Vec); 28*f2fddbb2SStefano Zampini PetscErrorCode apply_pc(PC,Vec,Vec); 29c4762a1bSJed Brown 30c4762a1bSJed Brown int main(int Argc,char **Args) 31c4762a1bSJed Brown { 32c4762a1bSJed Brown PetscInt x_mesh = 15,levels = 3,cycles = 1,use_jacobi = 0; 33c4762a1bSJed Brown PetscInt i,smooths = 1,*N,its; 34c4762a1bSJed Brown PetscErrorCode ierr; 35c4762a1bSJed Brown PCMGType am = PC_MG_MULTIPLICATIVE; 36c4762a1bSJed Brown Mat cmat,mat[20],fmat; 37c4762a1bSJed Brown KSP cksp,ksp[20],kspmg; 38c4762a1bSJed Brown PetscReal e[3]; /* l_2 error,max error, residual */ 39c4762a1bSJed Brown const char *shellname; 40c4762a1bSJed Brown Vec x,solution,X[20],R[20],B[20]; 41c4762a1bSJed Brown PC pcmg,pc; 42c4762a1bSJed Brown PetscBool flg; 43c4762a1bSJed Brown 44c4762a1bSJed Brown ierr = PetscInitialize(&Argc,&Args,(char*)0,help);if (ierr) return ierr; 45c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-x",&x_mesh,NULL);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-l",&levels,NULL);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-c",&cycles,NULL);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-smooths",&smooths,NULL);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-a",&flg);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown if (flg) am = PC_MG_ADDITIVE; 52c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-f",&flg);CHKERRQ(ierr); 53c4762a1bSJed Brown if (flg) am = PC_MG_FULL; 54c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-j",&flg);CHKERRQ(ierr); 55c4762a1bSJed Brown if (flg) use_jacobi = 1; 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscMalloc1(levels,&N);CHKERRQ(ierr); 58c4762a1bSJed Brown N[0] = x_mesh; 59c4762a1bSJed Brown for (i=1; i<levels; i++) { 60c4762a1bSJed Brown N[i] = N[i-1]/2; 61c4762a1bSJed Brown if (N[i] < 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Too many levels or N is not large enough"); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = Create1dLaplacian(N[levels-1],&cmat);CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown ierr = KSPCreate(PETSC_COMM_WORLD,&kspmg);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = KSPGetPC(kspmg,&pcmg);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = KSPSetFromOptions(kspmg);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = PCSetType(pcmg,PCMG);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PCMGSetLevels(pcmg,levels,NULL);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = PCMGSetType(pcmg,am);CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = PCMGGetCoarseSolve(pcmg,&cksp);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = KSPSetOperators(cksp,cmat,cmat);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = KSPGetPC(cksp,&pc);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = KSPSetType(cksp,KSPPREONLY);CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* zero is finest level */ 80c4762a1bSJed Brown for (i=0; i<levels-1; i++) { 81*f2fddbb2SStefano Zampini Mat dummy; 82*f2fddbb2SStefano Zampini 83*f2fddbb2SStefano Zampini ierr = PCMGSetResidual(pcmg,levels - 1 - i,residual,NULL);CHKERRQ(ierr); 84*f2fddbb2SStefano Zampini ierr = MatCreateShell(PETSC_COMM_WORLD,N[i+1],N[i],N[i+1],N[i],NULL,&mat[i]);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = MatShellSetOperation(mat[i],MATOP_MULT,(void (*)(void))restrct);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = MatShellSetOperation(mat[i],MATOP_MULT_TRANSPOSE_ADD,(void (*)(void))interpolate);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = PCMGSetInterpolation(pcmg,levels - 1 - i,mat[i]);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = PCMGSetRestriction(pcmg,levels - 1 - i,mat[i]);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = PCMGSetCycleTypeOnLevel(pcmg,levels - 1 - i,(PCMGCycleType)cycles);CHKERRQ(ierr); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* set smoother */ 92c4762a1bSJed Brown ierr = PCMGGetSmoother(pcmg,levels - 1 - i,&ksp[i]);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = KSPGetPC(ksp[i],&pc);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = PCSetType(pc,PCSHELL);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PCShellSetName(pc,"user_precond");CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = PCShellGetName(pc,&shellname);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"level=%D, PCShell name is %s\n",i,shellname);CHKERRQ(ierr); 98c4762a1bSJed Brown 99*f2fddbb2SStefano Zampini /* this is not used unless different options are passed to the solver */ 100*f2fddbb2SStefano Zampini ierr = MatCreateShell(PETSC_COMM_WORLD,N[i],N[i],N[i],N[i],NULL,&dummy);CHKERRQ(ierr); 101*f2fddbb2SStefano Zampini ierr = MatShellSetOperation(dummy,MATOP_MULT,(void (*)(void))amult);CHKERRQ(ierr); 102*f2fddbb2SStefano Zampini ierr = KSPSetOperators(ksp[i],dummy,dummy);CHKERRQ(ierr); 103*f2fddbb2SStefano Zampini ierr = MatDestroy(&dummy);CHKERRQ(ierr); 104*f2fddbb2SStefano Zampini 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown We override the matrix passed in by forcing it to use Richardson with 107c4762a1bSJed Brown a user provided application. This is non-standard and this practice 108c4762a1bSJed Brown should be avoided. 109c4762a1bSJed Brown */ 110*f2fddbb2SStefano Zampini ierr = PCShellSetApply(pc,apply_pc);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = PCShellSetApplyRichardson(pc,gauss_seidel);CHKERRQ(ierr); 112c4762a1bSJed Brown if (use_jacobi) { 113c4762a1bSJed Brown ierr = PCShellSetApplyRichardson(pc,jacobi_smoother);CHKERRQ(ierr); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown ierr = KSPSetType(ksp[i],KSPRICHARDSON);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = KSPSetInitialGuessNonzero(ksp[i],PETSC_TRUE);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = KSPSetTolerances(ksp[i],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,smooths);CHKERRQ(ierr); 118c4762a1bSJed Brown 119c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N[i],&x);CHKERRQ(ierr); 120c4762a1bSJed Brown 121c4762a1bSJed Brown X[levels - 1 - i] = x; 122c4762a1bSJed Brown if (i > 0) { 123c4762a1bSJed Brown ierr = PCMGSetX(pcmg,levels - 1 - i,x);CHKERRQ(ierr); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N[i],&x);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown B[levels -1 - i] = x; 128c4762a1bSJed Brown if (i > 0) { 129c4762a1bSJed Brown ierr = PCMGSetRhs(pcmg,levels - 1 - i,x);CHKERRQ(ierr); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N[i],&x);CHKERRQ(ierr); 132c4762a1bSJed Brown 133c4762a1bSJed Brown R[levels - 1 - i] = x; 134c4762a1bSJed Brown 135c4762a1bSJed Brown ierr = PCMGSetR(pcmg,levels - 1 - i,x);CHKERRQ(ierr); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown /* create coarse level vectors */ 138c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = PCMGSetX(pcmg,0,x);CHKERRQ(ierr); X[0] = x; 140c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = PCMGSetRhs(pcmg,0,x);CHKERRQ(ierr); B[0] = x; 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* create matrix multiply for finest level */ 144*f2fddbb2SStefano Zampini ierr = MatCreateShell(PETSC_COMM_WORLD,N[0],N[0],N[0],N[0],NULL,&fmat);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatShellSetOperation(fmat,MATOP_MULT,(void (*)(void))amult);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = KSPSetOperators(kspmg,fmat,fmat);CHKERRQ(ierr); 147c4762a1bSJed Brown 148c4762a1bSJed Brown ierr = CalculateSolution(N[0],&solution);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = CalculateRhs(B[levels-1]);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = VecSet(X[levels-1],0.0);CHKERRQ(ierr); 151c4762a1bSJed Brown 152c4762a1bSJed Brown ierr = residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = CalculateError(solution,X[levels-1],R[levels-1],e);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"l_2 error %g max error %g resi %g\n",(double)e[0],(double)e[1],(double)e[2]);CHKERRQ(ierr); 155c4762a1bSJed Brown 156c4762a1bSJed Brown ierr = KSPSolve(kspmg,B[levels-1],X[levels-1]);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = KSPGetIterationNumber(kspmg,&its);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = residual((Mat)0,B[levels-1],X[levels-1],R[levels-1]);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = CalculateError(solution,X[levels-1],R[levels-1],e);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"its %D l_2 error %g max error %g resi %g\n",its,(double)e[0],(double)e[1],(double)e[2]);CHKERRQ(ierr); 161c4762a1bSJed Brown 162c4762a1bSJed Brown ierr = PetscFree(N);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = VecDestroy(&solution);CHKERRQ(ierr); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* note we have to keep a list of all vectors allocated, this is 166c4762a1bSJed Brown not ideal, but putting it in MGDestroy is not so good either*/ 167c4762a1bSJed Brown for (i=0; i<levels; i++) { 168c4762a1bSJed Brown ierr = VecDestroy(&X[i]);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = VecDestroy(&B[i]);CHKERRQ(ierr); 170c4762a1bSJed Brown if (i) {ierr = VecDestroy(&R[i]);CHKERRQ(ierr);} 171c4762a1bSJed Brown } 172c4762a1bSJed Brown for (i=0; i<levels-1; i++) { 173c4762a1bSJed Brown ierr = MatDestroy(&mat[i]);CHKERRQ(ierr); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown ierr = MatDestroy(&cmat);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = MatDestroy(&fmat);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = KSPDestroy(&kspmg);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = PetscFinalize(); 179c4762a1bSJed Brown return ierr; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 183c4762a1bSJed Brown PetscErrorCode residual(Mat mat,Vec bb,Vec xx,Vec rr) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown PetscInt i,n1; 186c4762a1bSJed Brown PetscErrorCode ierr; 187c4762a1bSJed Brown PetscScalar *x,*r; 188c4762a1bSJed Brown const PetscScalar *b; 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscFunctionBegin; 191c4762a1bSJed Brown ierr = VecGetSize(bb,&n1);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 195c4762a1bSJed Brown n1--; 196c4762a1bSJed Brown r[0] = b[0] + x[1] - 2.0*x[0]; 197c4762a1bSJed Brown r[n1] = b[n1] + x[n1-1] - 2.0*x[n1]; 198c4762a1bSJed Brown for (i=1; i<n1; i++) r[i] = b[i] + x[i+1] + x[i-1] - 2.0*x[i]; 199c4762a1bSJed Brown ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 202c4762a1bSJed Brown PetscFunctionReturn(0); 203c4762a1bSJed Brown } 204*f2fddbb2SStefano Zampini 205c4762a1bSJed Brown PetscErrorCode amult(Mat mat,Vec xx,Vec yy) 206c4762a1bSJed Brown { 207c4762a1bSJed Brown PetscInt i,n1; 208c4762a1bSJed Brown PetscErrorCode ierr; 209c4762a1bSJed Brown PetscScalar *y; 210c4762a1bSJed Brown const PetscScalar *x; 211c4762a1bSJed Brown 212c4762a1bSJed Brown PetscFunctionBegin; 213c4762a1bSJed Brown ierr = VecGetSize(xx,&n1);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 216c4762a1bSJed Brown n1--; 217c4762a1bSJed Brown y[0] = -x[1] + 2.0*x[0]; 218c4762a1bSJed Brown y[n1] = -x[n1-1] + 2.0*x[n1]; 219c4762a1bSJed Brown for (i=1; i<n1; i++) y[i] = -x[i+1] - x[i-1] + 2.0*x[i]; 220c4762a1bSJed Brown ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 222c4762a1bSJed Brown PetscFunctionReturn(0); 223c4762a1bSJed Brown } 224*f2fddbb2SStefano Zampini 225c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 226*f2fddbb2SStefano Zampini PetscErrorCode apply_pc(PC pc,Vec bb,Vec xx) 227*f2fddbb2SStefano Zampini { 228*f2fddbb2SStefano Zampini PetscFunctionBegin; 229*f2fddbb2SStefano Zampini SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not implemented"); 230*f2fddbb2SStefano Zampini PetscFunctionReturn(0); 231*f2fddbb2SStefano Zampini } 232*f2fddbb2SStefano Zampini 233c4762a1bSJed Brown PetscErrorCode gauss_seidel(PC pc,Vec bb,Vec xx,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt m,PetscBool guesszero,PetscInt *its,PCRichardsonConvergedReason *reason) 234c4762a1bSJed Brown { 235c4762a1bSJed Brown PetscInt i,n1; 236c4762a1bSJed Brown PetscErrorCode ierr; 237c4762a1bSJed Brown PetscScalar *x; 238c4762a1bSJed Brown const PetscScalar *b; 239c4762a1bSJed Brown 240c4762a1bSJed Brown PetscFunctionBegin; 241c4762a1bSJed Brown *its = m; 242c4762a1bSJed Brown *reason = PCRICHARDSON_CONVERGED_ITS; 243c4762a1bSJed Brown ierr = VecGetSize(bb,&n1);CHKERRQ(ierr); n1--; 244c4762a1bSJed Brown ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 246c4762a1bSJed Brown while (m--) { 247c4762a1bSJed Brown x[0] = .5*(x[1] + b[0]); 248c4762a1bSJed Brown for (i=1; i<n1; i++) x[i] = .5*(x[i+1] + x[i-1] + b[i]); 249c4762a1bSJed Brown x[n1] = .5*(x[n1-1] + b[n1]); 250c4762a1bSJed Brown for (i=n1-1; i>0; i--) x[i] = .5*(x[i+1] + x[i-1] + b[i]); 251c4762a1bSJed Brown x[0] = .5*(x[1] + b[0]); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 255c4762a1bSJed Brown PetscFunctionReturn(0); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 258c4762a1bSJed Brown PetscErrorCode jacobi_smoother(PC pc,Vec bb,Vec xx,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt m,PetscBool guesszero,PetscInt *its,PCRichardsonConvergedReason *reason) 259c4762a1bSJed Brown { 260c4762a1bSJed Brown PetscInt i,n,n1; 261c4762a1bSJed Brown PetscErrorCode ierr; 262c4762a1bSJed Brown PetscScalar *r,*x; 263c4762a1bSJed Brown const PetscScalar *b; 264c4762a1bSJed Brown 265c4762a1bSJed Brown PetscFunctionBegin; 266c4762a1bSJed Brown *its = m; 267c4762a1bSJed Brown *reason = PCRICHARDSON_CONVERGED_ITS; 268c4762a1bSJed Brown ierr = VecGetSize(bb,&n);CHKERRQ(ierr); n1 = n - 1; 269c4762a1bSJed Brown ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = VecGetArray(w,&r);CHKERRQ(ierr); 272c4762a1bSJed Brown 273c4762a1bSJed Brown while (m--) { 274c4762a1bSJed Brown r[0] = .5*(x[1] + b[0]); 275c4762a1bSJed Brown for (i=1; i<n1; i++) r[i] = .5*(x[i+1] + x[i-1] + b[i]); 276c4762a1bSJed Brown r[n1] = .5*(x[n1-1] + b[n1]); 277c4762a1bSJed Brown for (i=0; i<n; i++) x[i] = (2.0*r[i] + x[i])/3.0; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 280c4762a1bSJed Brown ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 281c4762a1bSJed Brown ierr = VecRestoreArray(w,&r);CHKERRQ(ierr); 282c4762a1bSJed Brown PetscFunctionReturn(0); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown /* 285c4762a1bSJed Brown We know for this application that yy and zz are the same 286c4762a1bSJed Brown */ 287c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 288c4762a1bSJed Brown PetscErrorCode interpolate(Mat mat,Vec xx,Vec yy,Vec zz) 289c4762a1bSJed Brown { 290c4762a1bSJed Brown PetscInt i,n,N,i2; 291c4762a1bSJed Brown PetscErrorCode ierr; 292c4762a1bSJed Brown PetscScalar *y; 293c4762a1bSJed Brown const PetscScalar *x; 294c4762a1bSJed Brown 295c4762a1bSJed Brown PetscFunctionBegin; 296c4762a1bSJed Brown ierr = VecGetSize(yy,&N);CHKERRQ(ierr); 297c4762a1bSJed Brown ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 298c4762a1bSJed Brown ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 299c4762a1bSJed Brown n = N/2; 300c4762a1bSJed Brown for (i=0; i<n; i++) { 301c4762a1bSJed Brown i2 = 2*i; 302c4762a1bSJed Brown y[i2] += .5*x[i]; 303c4762a1bSJed Brown y[i2+1] += x[i]; 304c4762a1bSJed Brown y[i2+2] += .5*x[i]; 305c4762a1bSJed Brown } 306c4762a1bSJed Brown ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 308c4762a1bSJed Brown PetscFunctionReturn(0); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 311c4762a1bSJed Brown PetscErrorCode restrct(Mat mat,Vec rr,Vec bb) 312c4762a1bSJed Brown { 313c4762a1bSJed Brown PetscInt i,n,N,i2; 314c4762a1bSJed Brown PetscErrorCode ierr; 315c4762a1bSJed Brown PetscScalar *b; 316c4762a1bSJed Brown const PetscScalar *r; 317c4762a1bSJed Brown 318c4762a1bSJed Brown PetscFunctionBegin; 319c4762a1bSJed Brown ierr = VecGetSize(rr,&N);CHKERRQ(ierr); 320c4762a1bSJed Brown ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr); 321c4762a1bSJed Brown ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 322c4762a1bSJed Brown n = N/2; 323c4762a1bSJed Brown 324c4762a1bSJed Brown for (i=0; i<n; i++) { 325c4762a1bSJed Brown i2 = 2*i; 326c4762a1bSJed Brown b[i] = (r[i2] + 2.0*r[i2+1] + r[i2+2]); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr); 329c4762a1bSJed Brown ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 330c4762a1bSJed Brown PetscFunctionReturn(0); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 333c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt n,Mat *mat) 334c4762a1bSJed Brown { 335c4762a1bSJed Brown PetscScalar mone = -1.0,two = 2.0; 336c4762a1bSJed Brown PetscInt i,idx; 337c4762a1bSJed Brown PetscErrorCode ierr; 338c4762a1bSJed Brown 339c4762a1bSJed Brown PetscFunctionBegin; 340c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,mat);CHKERRQ(ierr); 341c4762a1bSJed Brown 342c4762a1bSJed Brown idx = n-1; 343c4762a1bSJed Brown ierr = MatSetValues(*mat,1,&idx,1,&idx,&two,INSERT_VALUES);CHKERRQ(ierr); 344c4762a1bSJed Brown for (i=0; i<n-1; i++) { 345c4762a1bSJed Brown ierr = MatSetValues(*mat,1,&i,1,&i,&two,INSERT_VALUES);CHKERRQ(ierr); 346c4762a1bSJed Brown idx = i+1; 347c4762a1bSJed Brown ierr = MatSetValues(*mat,1,&idx,1,&i,&mone,INSERT_VALUES);CHKERRQ(ierr); 348c4762a1bSJed Brown ierr = MatSetValues(*mat,1,&i,1,&idx,&mone,INSERT_VALUES);CHKERRQ(ierr); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351c4762a1bSJed Brown ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352c4762a1bSJed Brown PetscFunctionReturn(0); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 355c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec u) 356c4762a1bSJed Brown { 357c4762a1bSJed Brown PetscErrorCode ierr; 358c4762a1bSJed Brown PetscInt i,n; 359c4762a1bSJed Brown PetscReal h,x = 0.0; 360c4762a1bSJed Brown PetscScalar uu; 361c4762a1bSJed Brown 362c4762a1bSJed Brown PetscFunctionBegin; 363c4762a1bSJed Brown ierr = VecGetSize(u,&n);CHKERRQ(ierr); 364c4762a1bSJed Brown h = 1.0/((PetscReal)(n+1)); 365c4762a1bSJed Brown for (i=0; i<n; i++) { 366c4762a1bSJed Brown x += h; uu = 2.0*h*h; 367c4762a1bSJed Brown ierr = VecSetValues(u,1,&i,&uu,INSERT_VALUES);CHKERRQ(ierr); 368c4762a1bSJed Brown } 369c4762a1bSJed Brown PetscFunctionReturn(0); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 372c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt n,Vec *solution) 373c4762a1bSJed Brown { 374c4762a1bSJed Brown PetscErrorCode ierr; 375c4762a1bSJed Brown PetscInt i; 376c4762a1bSJed Brown PetscReal h,x = 0.0; 377c4762a1bSJed Brown PetscScalar uu; 378c4762a1bSJed Brown 379c4762a1bSJed Brown PetscFunctionBegin; 380c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,solution);CHKERRQ(ierr); 381c4762a1bSJed Brown h = 1.0/((PetscReal)(n+1)); 382c4762a1bSJed Brown for (i=0; i<n; i++) { 383c4762a1bSJed Brown x += h; uu = x*(1.-x); 384c4762a1bSJed Brown ierr = VecSetValues(*solution,1,&i,&uu,INSERT_VALUES);CHKERRQ(ierr); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown PetscFunctionReturn(0); 387c4762a1bSJed Brown } 388c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 389c4762a1bSJed Brown PetscErrorCode CalculateError(Vec solution,Vec u,Vec r,PetscReal *e) 390c4762a1bSJed Brown { 391c4762a1bSJed Brown PetscErrorCode ierr; 392c4762a1bSJed Brown 393c4762a1bSJed Brown PetscFunctionBegin; 394c4762a1bSJed Brown ierr = VecNorm(r,NORM_2,e+2);CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = VecWAXPY(r,-1.0,u,solution);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = VecNorm(r,NORM_2,e);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = VecNorm(r,NORM_1,e+1);CHKERRQ(ierr); 398c4762a1bSJed Brown PetscFunctionReturn(0); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown 401c4762a1bSJed Brown /*TEST 402c4762a1bSJed Brown 403c4762a1bSJed Brown test: 404c4762a1bSJed Brown 405c4762a1bSJed Brown TEST*/ 406