Actual source code: ex66.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: /*   DMDA/KSP solving a system of linear equations.
  2:      Poisson equation in 2D:

  4:      div(grad p) = f,  0 < x,y < 1
  5:      with
  6:        forcing function f = -cos(m*pi*x)*cos(n*pi*y),
  7:        Periodic boundary conditions
  8:          p(x=0) = p(x=1)
  9:        Neuman boundary conditions
 10:          dp/dy = 0 for y = 0, y = 1.

 12:        This example uses the DM_BOUNDARY_MIRROR to implement the Neumann boundary conditions, see the manual pages for DMBoundaryType

 14:        Compare to ex50.c
 15: */

 17: static char help[] = "Solves 2D Poisson equation using multigrid.\n\n";

 19:  #include <petscdm.h>
 20:  #include <petscdmda.h>
 21:  #include <petscksp.h>
 22:  #include <petscsys.h>
 23:  #include <petscvec.h>

 25: extern PetscErrorCode ComputeJacobian(KSP,Mat,Mat,void*);
 26: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);

 28: typedef struct {
 29:   PetscScalar uu, tt;
 30: } UserContext;

 32: int main(int argc,char **argv)
 33: {
 34:   KSP            ksp;
 35:   DM             da;
 36:   UserContext    user;

 39:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 40:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 41:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_MIRROR,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 42:   DMSetFromOptions(da);
 43:   DMSetUp(da);
 44:   KSPSetDM(ksp,(DM)da);
 45:   DMSetApplicationContext(da,&user);

 47:   user.uu     = 1.0;
 48:   user.tt     = 1.0;

 50:   KSPSetComputeRHS(ksp,ComputeRHS,&user);
 51:   KSPSetComputeOperators(ksp,ComputeJacobian,&user);
 52:   KSPSetFromOptions(ksp);
 53:   KSPSolve(ksp,NULL,NULL);

 55:   DMDestroy(&da);
 56:   KSPDestroy(&ksp);
 57:   PetscFinalize();
 58:   return ierr;
 59: }

 61: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
 62: {
 63:   UserContext    *user = (UserContext*)ctx;
 65:   PetscInt       i,j,M,N,xm,ym,xs,ys;
 66:   PetscScalar    Hx,Hy,pi,uu,tt;
 67:   PetscScalar    **array;
 68:   DM             da;
 69:   MatNullSpace   nullspace;

 72:   KSPGetDM(ksp,&da);
 73:   DMDAGetInfo(da, 0, &M, &N, 0,0,0,0,0,0,0,0,0,0);
 74:   uu   = user->uu; tt = user->tt;
 75:   pi   = 4*PetscAtanReal(1.0);
 76:   Hx   = 1.0/(PetscReal)(M);
 77:   Hy   = 1.0/(PetscReal)(N);

 79:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0); /* Fine grid */
 80:   DMDAVecGetArray(da, b, &array);
 81:   for (j=ys; j<ys+ym; j++) {
 82:     for (i=xs; i<xs+xm; i++) {
 83:       array[j][i] = -PetscCosScalar(uu*pi*((PetscReal)i+0.5)*Hx)*PetscCosScalar(tt*pi*((PetscReal)j+0.5)*Hy)*Hx*Hy;
 84:     }
 85:   }
 86:   DMDAVecRestoreArray(da, b, &array);
 87:   VecAssemblyBegin(b);
 88:   VecAssemblyEnd(b);

 90:   /* force right hand side to be consistent for singular matrix */
 91:   /* note this is really a hack, normally the model would provide you with a consistent right handside */
 92:   MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
 93:   MatNullSpaceRemove(nullspace,b);
 94:   MatNullSpaceDestroy(&nullspace);
 95:   return(0);
 96: }

 98: PetscErrorCode ComputeJacobian(KSP ksp,Mat J, Mat jac,void *ctx)
 99: {
101:   PetscInt       i, j, M, N, xm, ym, xs, ys;
102:   PetscScalar    v[5], Hx, Hy, HydHx, HxdHy;
103:   MatStencil     row, col[5];
104:   DM             da;
105:   MatNullSpace   nullspace;

108:   KSPGetDM(ksp,&da);
109:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
110:   Hx    = 1.0 / (PetscReal)(M);
111:   Hy    = 1.0 / (PetscReal)(N);
112:   HxdHy = Hx/Hy;
113:   HydHx = Hy/Hx;
114:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
115:   for (j=ys; j<ys+ym; j++) {
116:     for (i=xs; i<xs+xm; i++) {
117:       row.i = i; row.j = j;
118:       v[0] = -HxdHy;              col[0].i = i;   col[0].j = j-1;
119:       v[1] = -HydHx;              col[1].i = i-1; col[1].j = j;
120:       v[2] = 2.0*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
121:       v[3] = -HydHx;              col[3].i = i+1; col[3].j = j;
122:       v[4] = -HxdHy;              col[4].i = i;   col[4].j = j+1;
123:       MatSetValuesStencil(jac,1,&row,5,col,v,ADD_VALUES);
124:     }
125:   }
126:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
127:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

129:   MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
130:   MatSetNullSpace(J,nullspace);
131:   MatNullSpaceDestroy(&nullspace);
132:   return(0);
133: }



137: /*TEST

139:    build:
140:       requires: !complex !single

142:    test:
143:       args: -pc_type mg -pc_mg_type full -ksp_monitor_short -da_refine 3 -mg_coarse_pc_type svd -ksp_view

145:    test:
146:       suffix: 2
147:       nsize: 4
148:       args: -pc_type mg -pc_mg_type full -ksp_monitor_short -da_refine 3 -mg_coarse_pc_type redundant -mg_coarse_redundant_pc_type svd -ksp_view

150: TEST*/