Actual source code: ex28.c

petsc-3.10.5 2019-03-28
Report Typos and Errors


  3: static char help[] = "Solves 1D wave equation using multigrid.\n\n";

  5:  #include <petscdm.h>
  6:  #include <petscdmda.h>
  7:  #include <petscksp.h>


 10: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
 11: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
 12: extern PetscErrorCode ComputeInitialSolution(DM,Vec);

 14: int main(int argc,char **argv)
 15: {
 17:   PetscInt       i;
 18:   KSP            ksp;
 19:   DM             da;
 20:   Vec            x;

 22:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 23:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 24:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,3,2,1,0,&da);
 25:   DMSetFromOptions(da);
 26:   DMSetUp(da);
 27:   KSPSetDM(ksp,da);
 28:   KSPSetComputeRHS(ksp,ComputeRHS,NULL);
 29:   KSPSetComputeOperators(ksp,ComputeMatrix,NULL);

 31:   KSPSetFromOptions(ksp);
 32:   DMCreateGlobalVector(da,&x);
 33:   ComputeInitialSolution(da,x);
 34:   DMSetApplicationContext(da,x);
 35:   KSPSetUp(ksp);
 36:   VecView(x,PETSC_VIEWER_DRAW_WORLD);
 37:   for (i=0; i<10; i++) {
 38:     KSPSolve(ksp,NULL,x);
 39:     VecView(x,PETSC_VIEWER_DRAW_WORLD);
 40:   }
 41:   VecDestroy(&x);
 42:   KSPDestroy(&ksp);
 43:   DMDestroy(&da);
 44:   PetscFinalize();
 45:   return ierr;
 46: }

 48: PetscErrorCode ComputeInitialSolution(DM da,Vec x)
 49: {
 51:   PetscInt       mx,col[2],xs,xm,i;
 52:   PetscScalar    Hx,val[2];

 55:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
 56:   Hx   = 2.0*PETSC_PI / (PetscReal)(mx);
 57:   DMDAGetCorners(da,&xs,0,0,&xm,0,0);

 59:   for (i=xs; i<xs+xm; i++) {
 60:     col[0] = 2*i; col[1] = 2*i + 1;
 61:     val[0] = val[1] = PetscSinScalar(((PetscScalar)i)*Hx);
 62:     VecSetValues(x,2,col,val,INSERT_VALUES);
 63:   }
 64:   VecAssemblyBegin(x);
 65:   VecAssemblyEnd(x);
 66:   return(0);
 67: }

 69: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
 70: {
 72:   PetscInt       mx;
 73:   PetscScalar    h;
 74:   Vec            x;
 75:   DM             da;

 78:   KSPGetDM(ksp,&da);
 79:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
 80:   DMGetApplicationContext(da,&x);
 81:   h    = 2.0*PETSC_PI/((mx));
 82:   VecCopy(x,b);
 83:   VecScale(b,h);
 84:   return(0);
 85: }

 87: PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
 88: {
 90:   PetscInt       i,mx,xm,xs;
 91:   PetscScalar    v[7],Hx;
 92:   MatStencil     row,col[7];
 93:   PetscScalar    lambda;
 94:   DM             da;

 97:   KSPGetDM(ksp,&da);
 98:   PetscMemzero(col,7*sizeof(MatStencil));
 99:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
100:   Hx     = 2.0*PETSC_PI / (PetscReal)(mx);
101:   DMDAGetCorners(da,&xs,0,0,&xm,0,0);
102:   lambda = 2.0*Hx;
103:   for (i=xs; i<xs+xm; i++) {
104:     row.i = i; row.j = 0; row.k = 0; row.c = 0;
105:     v[0]  = Hx;     col[0].i = i;   col[0].c = 0;
106:     v[1]  = lambda; col[1].i = i-1;   col[1].c = 1;
107:     v[2]  = -lambda;col[2].i = i+1; col[2].c = 1;
108:     MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);

110:     row.i = i; row.j = 0; row.k = 0; row.c = 1;
111:     v[0]  = lambda; col[0].i = i-1;   col[0].c = 0;
112:     v[1]  = Hx;     col[1].i = i;   col[1].c = 1;
113:     v[2]  = -lambda;col[2].i = i+1; col[2].c = 0;
114:     MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
115:   }
116:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
117:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
118:   MatView(jac,PETSC_VIEWER_BINARY_(PETSC_COMM_SELF));
119:   return(0);
120: }


123: /*TEST

125:    test:
126:       args: -ksp_monitor_short -pc_type mg -pc_mg_type full -ksp_type fgmres -da_refine 2 -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 1 -mg_levels_pc_type ilu

128: TEST*/