Actual source code: ex25.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: /*
  3:  Partial differential equation

  5:    d  (1 + e*sine(2*pi*k*x)) d u = 1, 0 < x < 1,
  6:    --                        ---
  7:    dx                        dx
  8: with boundary conditions

 10:    u = 0 for x = 0, x = 1

 12:    This uses multigrid to solve the linear system

 14: */

 16: static char help[] = "Solves 1D variable coefficient Laplacian using multigrid.\n\n";

 18: #include <petscdm.h>
 19: #include <petscdmda.h>
 20: #include <petscksp.h>

 22: static PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
 23: static PetscErrorCode ComputeRHS(KSP,Vec,void*);

 25: typedef struct {
 26:   PetscInt    k;
 27:   PetscScalar e;
 28: } AppCtx;

 32: int main(int argc,char **argv)
 33: {
 35:   KSP            ksp;
 36:   DM             da;
 37:   AppCtx         user;
 38:   Mat            A;
 39:   Vec            b,b2;
 40:   Vec            x;
 41:   PetscReal      nrm;

 43:   PetscInitialize(&argc,&argv,(char*)0,help);

 45:   user.k = 1;
 46:   user.e = .99;
 47:   PetscOptionsGetInt(NULL,0,"-k",&user.k,0);
 48:   PetscOptionsGetScalar(NULL,0,"-e",&user.e,0);

 50:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 51:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,-128,1,1,0,&da);
 52:   KSPSetDM(ksp,da);
 53:   KSPSetComputeRHS(ksp,ComputeRHS,&user);
 54:   KSPSetComputeOperators(ksp,ComputeMatrix,&user);
 55:   KSPSetFromOptions(ksp);
 56:   KSPSolve(ksp,NULL,NULL);

 58:   KSPGetOperators(ksp,&A,NULL);
 59:   KSPGetSolution(ksp,&x);
 60:   KSPGetRhs(ksp,&b);
 61:   VecDuplicate(b,&b2);
 62:   MatMult(A,x,b2);
 63:   VecAXPY(b2,-1.0,b);
 64:   VecNorm(b2,NORM_MAX,&nrm);
 65:   PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)nrm);

 67:   VecDestroy(&b2);
 68:   KSPDestroy(&ksp);
 69:   DMDestroy(&da);
 70:   PetscFinalize();

 72:   return 0;
 73: }

 77: static PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
 78: {
 80:   PetscInt       mx,idx[2];
 81:   PetscScalar    h,v[2];
 82:   DM             da;

 85:   KSPGetDM(ksp,&da);
 86:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
 87:   h      = 1.0/((mx-1));
 88:   VecSet(b,h);
 89:   idx[0] = 0; idx[1] = mx -1;
 90:   v[0]   = v[1] = 0.0;
 91:   VecSetValues(b,2,idx,v,INSERT_VALUES);
 92:   VecAssemblyBegin(b);
 93:   VecAssemblyEnd(b);
 94:   return(0);
 95: }

 99: static PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
100: {
101:   AppCtx         *user = (AppCtx*)ctx;
103:   PetscInt       i,mx,xm,xs;
104:   PetscScalar    v[3],h,xlow,xhigh;
105:   MatStencil     row,col[3];
106:   DM             da;

109:   KSPGetDM(ksp,&da);
110:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
111:   DMDAGetCorners(da,&xs,0,0,&xm,0,0);
112:   h    = 1.0/(mx-1);

114:   for (i=xs; i<xs+xm; i++) {
115:     row.i = i;
116:     if (i==0 || i==mx-1) {
117:       v[0] = 2.0/h;
118:       MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
119:     } else {
120:       xlow  = h*(PetscReal)i - .5*h;
121:       xhigh = xlow + h;
122:       v[0]  = (-1.0 - user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xlow))/h;col[0].i = i-1;
123:       v[1]  = (2.0 + user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xlow) + user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xhigh))/h;col[1].i = row.i;
124:       v[2]  = (-1.0 - user->e*PetscSinScalar(2.0*PETSC_PI*user->k*xhigh))/h;col[2].i = i+1;
125:       MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
126:     }
127:   }
128:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
129:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
130:   return(0);
131: }