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); 28f2fddbb2SStefano 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; 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-x",&x_mesh,NULL)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-l",&levels,NULL)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-c",&cycles,NULL)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-smooths",&smooths,NULL)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-a",&flg)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown if (flg) am = PC_MG_ADDITIVE; 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-f",&flg)); 53c4762a1bSJed Brown if (flg) am = PC_MG_FULL; 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-j",&flg)); 55c4762a1bSJed Brown if (flg) use_jacobi = 1; 56c4762a1bSJed Brown 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(levels,&N)); 58c4762a1bSJed Brown N[0] = x_mesh; 59c4762a1bSJed Brown for (i=1; i<levels; i++) { 60c4762a1bSJed Brown N[i] = N[i-1]/2; 612c71b3e2SJacob Faibussowitsch PetscCheckFalse(N[i] < 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"Too many levels or N is not large enough"); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(Create1dLaplacian(N[levels-1],&cmat)); 65c4762a1bSJed Brown 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&kspmg)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(kspmg,&pcmg)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(kspmg)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pcmg,PCMG)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetLevels(pcmg,levels,NULL)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetType(pcmg,am)); 72c4762a1bSJed Brown 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetCoarseSolve(pcmg,&cksp)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(cksp,cmat,cmat)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(cksp,&pc)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCLU)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(cksp,KSPPREONLY)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* zero is finest level */ 80c4762a1bSJed Brown for (i=0; i<levels-1; i++) { 81f2fddbb2SStefano Zampini Mat dummy; 82f2fddbb2SStefano Zampini 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetResidual(pcmg,levels - 1 - i,residual,NULL)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,N[i+1],N[i],N[i+1],N[i],NULL,&mat[i])); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(mat[i],MATOP_MULT,(void (*)(void))restrct)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(mat[i],MATOP_MULT_TRANSPOSE_ADD,(void (*)(void))interpolate)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetInterpolation(pcmg,levels - 1 - i,mat[i])); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetRestriction(pcmg,levels - 1 - i,mat[i])); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetCycleTypeOnLevel(pcmg,levels - 1 - i,(PCMGCycleType)cycles)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* set smoother */ 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetSmoother(pcmg,levels - 1 - i,&ksp[i])); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp[i],&pc)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCSHELL)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCShellSetName(pc,"user_precond")); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCShellGetName(pc,&shellname)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"level=%D, PCShell name is %s\n",i,shellname)); 98c4762a1bSJed Brown 99f2fddbb2SStefano Zampini /* this is not used unless different options are passed to the solver */ 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,N[i],N[i],N[i],N[i],NULL,&dummy)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(dummy,MATOP_MULT,(void (*)(void))amult)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ksp[i],dummy,dummy)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dummy)); 104f2fddbb2SStefano 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*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCShellSetApply(pc,apply_pc)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCShellSetApplyRichardson(pc,gauss_seidel)); 112c4762a1bSJed Brown if (use_jacobi) { 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCShellSetApplyRichardson(pc,jacobi_smoother)); 114c4762a1bSJed Brown } 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(ksp[i],KSPRICHARDSON)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetInitialGuessNonzero(ksp[i],PETSC_TRUE)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetTolerances(ksp[i],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,smooths)); 118c4762a1bSJed Brown 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[i],&x)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown X[levels - 1 - i] = x; 122c4762a1bSJed Brown if (i > 0) { 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetX(pcmg,levels - 1 - i,x)); 124c4762a1bSJed Brown } 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[i],&x)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown B[levels -1 - i] = x; 128c4762a1bSJed Brown if (i > 0) { 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetRhs(pcmg,levels - 1 - i,x)); 130c4762a1bSJed Brown } 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[i],&x)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown R[levels - 1 - i] = x; 134c4762a1bSJed Brown 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetR(pcmg,levels - 1 - i,x)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown /* create coarse level vectors */ 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetX(pcmg,0,x)); X[0] = x; 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N[levels-1],&x)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetRhs(pcmg,0,x)); B[0] = x; 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* create matrix multiply for finest level */ 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,N[0],N[0],N[0],N[0],NULL,&fmat)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(fmat,MATOP_MULT,(void (*)(void))amult)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(kspmg,fmat,fmat)); 147c4762a1bSJed Brown 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(CalculateSolution(N[0],&solution)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(CalculateRhs(B[levels-1])); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X[levels-1],0.0)); 151c4762a1bSJed Brown 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(residual((Mat)0,B[levels-1],X[levels-1],R[levels-1])); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(CalculateError(solution,X[levels-1],R[levels-1],e)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"l_2 error %g max error %g resi %g\n",(double)e[0],(double)e[1],(double)e[2])); 155c4762a1bSJed Brown 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(kspmg,B[levels-1],X[levels-1])); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetIterationNumber(kspmg,&its)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(residual((Mat)0,B[levels-1],X[levels-1],R[levels-1])); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(CalculateError(solution,X[levels-1],R[levels-1],e)); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(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])); 161c4762a1bSJed Brown 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(N)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&solution)); 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++) { 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X[i])); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&B[i])); 170*5f80ce2aSJacob Faibussowitsch if (i) CHKERRQ(VecDestroy(&R[i])); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown for (i=0; i<levels-1; i++) { 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat[i])); 174c4762a1bSJed Brown } 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&cmat)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&fmat)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&kspmg)); 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 PetscScalar *x,*r; 187c4762a1bSJed Brown const PetscScalar *b; 188c4762a1bSJed Brown 189c4762a1bSJed Brown PetscFunctionBegin; 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(bb,&n1)); 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(bb,&b)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xx,&x)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(rr,&r)); 194c4762a1bSJed Brown n1--; 195c4762a1bSJed Brown r[0] = b[0] + x[1] - 2.0*x[0]; 196c4762a1bSJed Brown r[n1] = b[n1] + x[n1-1] - 2.0*x[n1]; 197c4762a1bSJed Brown for (i=1; i<n1; i++) r[i] = b[i] + x[i+1] + x[i-1] - 2.0*x[i]; 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(bb,&b)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xx,&x)); 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(rr,&r)); 201c4762a1bSJed Brown PetscFunctionReturn(0); 202c4762a1bSJed Brown } 203f2fddbb2SStefano Zampini 204c4762a1bSJed Brown PetscErrorCode amult(Mat mat,Vec xx,Vec yy) 205c4762a1bSJed Brown { 206c4762a1bSJed Brown PetscInt i,n1; 207c4762a1bSJed Brown PetscScalar *y; 208c4762a1bSJed Brown const PetscScalar *x; 209c4762a1bSJed Brown 210c4762a1bSJed Brown PetscFunctionBegin; 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(xx,&n1)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(xx,&x)); 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(yy,&y)); 214c4762a1bSJed Brown n1--; 215c4762a1bSJed Brown y[0] = -x[1] + 2.0*x[0]; 216c4762a1bSJed Brown y[n1] = -x[n1-1] + 2.0*x[n1]; 217c4762a1bSJed Brown for (i=1; i<n1; i++) y[i] = -x[i+1] - x[i-1] + 2.0*x[i]; 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(xx,&x)); 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(yy,&y)); 220c4762a1bSJed Brown PetscFunctionReturn(0); 221c4762a1bSJed Brown } 222f2fddbb2SStefano Zampini 223c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 224f2fddbb2SStefano Zampini PetscErrorCode apply_pc(PC pc,Vec bb,Vec xx) 225f2fddbb2SStefano Zampini { 226f2fddbb2SStefano Zampini PetscFunctionBegin; 227f2fddbb2SStefano Zampini SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not implemented"); 228f2fddbb2SStefano Zampini } 229f2fddbb2SStefano Zampini 230c4762a1bSJed 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) 231c4762a1bSJed Brown { 232c4762a1bSJed Brown PetscInt i,n1; 233c4762a1bSJed Brown PetscScalar *x; 234c4762a1bSJed Brown const PetscScalar *b; 235c4762a1bSJed Brown 236c4762a1bSJed Brown PetscFunctionBegin; 237c4762a1bSJed Brown *its = m; 238c4762a1bSJed Brown *reason = PCRICHARDSON_CONVERGED_ITS; 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(bb,&n1)); n1--; 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(bb,&b)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xx,&x)); 242c4762a1bSJed Brown while (m--) { 243c4762a1bSJed Brown x[0] = .5*(x[1] + b[0]); 244c4762a1bSJed Brown for (i=1; i<n1; i++) x[i] = .5*(x[i+1] + x[i-1] + b[i]); 245c4762a1bSJed Brown x[n1] = .5*(x[n1-1] + b[n1]); 246c4762a1bSJed Brown for (i=n1-1; i>0; i--) x[i] = .5*(x[i+1] + x[i-1] + b[i]); 247c4762a1bSJed Brown x[0] = .5*(x[1] + b[0]); 248c4762a1bSJed Brown } 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(bb,&b)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xx,&x)); 251c4762a1bSJed Brown PetscFunctionReturn(0); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 254c4762a1bSJed 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) 255c4762a1bSJed Brown { 256c4762a1bSJed Brown PetscInt i,n,n1; 257c4762a1bSJed Brown PetscScalar *r,*x; 258c4762a1bSJed Brown const PetscScalar *b; 259c4762a1bSJed Brown 260c4762a1bSJed Brown PetscFunctionBegin; 261c4762a1bSJed Brown *its = m; 262c4762a1bSJed Brown *reason = PCRICHARDSON_CONVERGED_ITS; 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(bb,&n)); n1 = n - 1; 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(bb,&b)); 265*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xx,&x)); 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(w,&r)); 267c4762a1bSJed Brown 268c4762a1bSJed Brown while (m--) { 269c4762a1bSJed Brown r[0] = .5*(x[1] + b[0]); 270c4762a1bSJed Brown for (i=1; i<n1; i++) r[i] = .5*(x[i+1] + x[i-1] + b[i]); 271c4762a1bSJed Brown r[n1] = .5*(x[n1-1] + b[n1]); 272c4762a1bSJed Brown for (i=0; i<n; i++) x[i] = (2.0*r[i] + x[i])/3.0; 273c4762a1bSJed Brown } 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(bb,&b)); 275*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xx,&x)); 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(w,&r)); 277c4762a1bSJed Brown PetscFunctionReturn(0); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown We know for this application that yy and zz are the same 281c4762a1bSJed Brown */ 282c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 283c4762a1bSJed Brown PetscErrorCode interpolate(Mat mat,Vec xx,Vec yy,Vec zz) 284c4762a1bSJed Brown { 285c4762a1bSJed Brown PetscInt i,n,N,i2; 286c4762a1bSJed Brown PetscScalar *y; 287c4762a1bSJed Brown const PetscScalar *x; 288c4762a1bSJed Brown 289c4762a1bSJed Brown PetscFunctionBegin; 290*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(yy,&N)); 291*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(xx,&x)); 292*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(yy,&y)); 293c4762a1bSJed Brown n = N/2; 294c4762a1bSJed Brown for (i=0; i<n; i++) { 295c4762a1bSJed Brown i2 = 2*i; 296c4762a1bSJed Brown y[i2] += .5*x[i]; 297c4762a1bSJed Brown y[i2+1] += x[i]; 298c4762a1bSJed Brown y[i2+2] += .5*x[i]; 299c4762a1bSJed Brown } 300*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(xx,&x)); 301*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(yy,&y)); 302c4762a1bSJed Brown PetscFunctionReturn(0); 303c4762a1bSJed Brown } 304c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 305c4762a1bSJed Brown PetscErrorCode restrct(Mat mat,Vec rr,Vec bb) 306c4762a1bSJed Brown { 307c4762a1bSJed Brown PetscInt i,n,N,i2; 308c4762a1bSJed Brown PetscScalar *b; 309c4762a1bSJed Brown const PetscScalar *r; 310c4762a1bSJed Brown 311c4762a1bSJed Brown PetscFunctionBegin; 312*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(rr,&N)); 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(rr,&r)); 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(bb,&b)); 315c4762a1bSJed Brown n = N/2; 316c4762a1bSJed Brown 317c4762a1bSJed Brown for (i=0; i<n; i++) { 318c4762a1bSJed Brown i2 = 2*i; 319c4762a1bSJed Brown b[i] = (r[i2] + 2.0*r[i2+1] + r[i2+2]); 320c4762a1bSJed Brown } 321*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(rr,&r)); 322*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(bb,&b)); 323c4762a1bSJed Brown PetscFunctionReturn(0); 324c4762a1bSJed Brown } 325c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 326c4762a1bSJed Brown PetscErrorCode Create1dLaplacian(PetscInt n,Mat *mat) 327c4762a1bSJed Brown { 328c4762a1bSJed Brown PetscScalar mone = -1.0,two = 2.0; 329c4762a1bSJed Brown PetscInt i,idx; 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscFunctionBegin; 332*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,mat)); 333c4762a1bSJed Brown 334c4762a1bSJed Brown idx = n-1; 335*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*mat,1,&idx,1,&idx,&two,INSERT_VALUES)); 336c4762a1bSJed Brown for (i=0; i<n-1; i++) { 337*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*mat,1,&i,1,&i,&two,INSERT_VALUES)); 338c4762a1bSJed Brown idx = i+1; 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*mat,1,&idx,1,&i,&mone,INSERT_VALUES)); 340*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(*mat,1,&i,1,&idx,&mone,INSERT_VALUES)); 341c4762a1bSJed Brown } 342*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY)); 343*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY)); 344c4762a1bSJed Brown PetscFunctionReturn(0); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 347c4762a1bSJed Brown PetscErrorCode CalculateRhs(Vec u) 348c4762a1bSJed Brown { 349c4762a1bSJed Brown PetscInt i,n; 350*5f80ce2aSJacob Faibussowitsch PetscReal h; 351c4762a1bSJed Brown PetscScalar uu; 352c4762a1bSJed Brown 353c4762a1bSJed Brown PetscFunctionBegin; 354*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(u,&n)); 355c4762a1bSJed Brown h = 1.0/((PetscReal)(n+1)); 356c4762a1bSJed Brown for (i=0; i<n; i++) { 357*5f80ce2aSJacob Faibussowitsch uu = 2.0*h*h; 358*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(u,1,&i,&uu,INSERT_VALUES)); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown PetscFunctionReturn(0); 361c4762a1bSJed Brown } 362c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 363c4762a1bSJed Brown PetscErrorCode CalculateSolution(PetscInt n,Vec *solution) 364c4762a1bSJed Brown { 365c4762a1bSJed Brown PetscInt i; 366c4762a1bSJed Brown PetscReal h,x = 0.0; 367c4762a1bSJed Brown PetscScalar uu; 368c4762a1bSJed Brown 369c4762a1bSJed Brown PetscFunctionBegin; 370*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,solution)); 371c4762a1bSJed Brown h = 1.0/((PetscReal)(n+1)); 372c4762a1bSJed Brown for (i=0; i<n; i++) { 373c4762a1bSJed Brown x += h; uu = x*(1.-x); 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(*solution,1,&i,&uu,INSERT_VALUES)); 375c4762a1bSJed Brown } 376c4762a1bSJed Brown PetscFunctionReturn(0); 377c4762a1bSJed Brown } 378c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 379c4762a1bSJed Brown PetscErrorCode CalculateError(Vec solution,Vec u,Vec r,PetscReal *e) 380c4762a1bSJed Brown { 381c4762a1bSJed Brown PetscFunctionBegin; 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r,NORM_2,e+2)); 383*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(r,-1.0,u,solution)); 384*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r,NORM_2,e)); 385*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r,NORM_1,e+1)); 386c4762a1bSJed Brown PetscFunctionReturn(0); 387c4762a1bSJed Brown } 388c4762a1bSJed Brown 389c4762a1bSJed Brown /*TEST 390c4762a1bSJed Brown 391c4762a1bSJed Brown test: 392c4762a1bSJed Brown 393c4762a1bSJed Brown TEST*/ 394