Actual source code: ex39.c

petsc-3.7.3 2016-07-24
Report Typos and Errors
  1: /*
  2: mpiexec -n 8 ./ex39 -ksp_type fbcgs -ksp_rtol 1.e-6 -sub_ksp_type bcgs -sub_ksp_rtol 1.e-3 -pc_type bjacobi -ksp_converged_reason -ksp_monitor -n1 32 -n2 32 -n3 32

  4:   Contributed by Jie Chen for testing flexible BiCGStab algorithm
  5: */

  7: static char help[] = "Solves the PDE (in 3D) - laplacian(u) + gamma x dot grad(u) + beta u = 1\n\
  8: with zero Dirichlet condition. The discretization is standard centered\n\
  9: difference. Input parameters include:\n\
 10:   -n1        : number of mesh points in 1st dimension (default 32)\n\
 11:   -n2        : number of mesh points in 2nd dimension (default 32)\n\
 12:   -n3        : number of mesh points in 3nd dimension (default 32)\n\
 13:   -h         : spacing between mesh points (default 1/n1)\n\
 14:   -gamma     : gamma (default 4/h)\n\
 15:   -beta      : beta (default 0.01/h^2)\n\n";

 17: #include <petscksp.h>
 20: int main(int argc,char **args)
 21: {
 22:   Vec            x,b,u;                 /* approx solution, RHS, working vector */
 23:   Mat            A;                     /* linear system matrix */
 24:   KSP            ksp;                   /* linear solver context */
 25:   PetscInt       n1, n2, n3;            /* parameters */
 26:   PetscReal      h, gamma, beta;        /* parameters */
 27:   PetscInt       i,j,k,Ii,J,Istart,Iend;
 29:   PetscScalar    v, co1, co2;

 31:   PetscInitialize(&argc,&args,(char*)0,help);
 32:   n1 = 32;
 33:   n2 = 32;
 34:   n3 = 32;
 35:   PetscOptionsGetInt(NULL,NULL,"-n1",&n1,NULL);
 36:   PetscOptionsGetInt(NULL,NULL,"-n2",&n2,NULL);
 37:   PetscOptionsGetInt(NULL,NULL,"-n3",&n3,NULL);

 39:   h     = 1.0/n1;
 40:   gamma = 4.0/h;
 41:   beta  = 0.01/(h*h);
 42:   PetscOptionsGetReal(NULL,NULL,"-h",&h,NULL);
 43:   PetscOptionsGetReal(NULL,NULL,"-gamma",&gamma,NULL);
 44:   PetscOptionsGetReal(NULL,NULL,"-beta",&beta,NULL);

 46:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 47:          Compute the matrix and set right-hand-side vector.
 48:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 49:   MatCreate(PETSC_COMM_WORLD,&A);
 50:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n1*n2*n3,n1*n2*n3);
 51:   MatSetFromOptions(A);
 52:   MatMPIAIJSetPreallocation(A,7,NULL,7,NULL);
 53:   MatSeqAIJSetPreallocation(A,7,NULL);
 54:   MatSetUp(A);
 55:   MatGetOwnershipRange(A,&Istart,&Iend);

 57:   /*
 58:      Set matrix elements for the 3-D, seven-point stencil in parallel.
 59:       - Each processor needs to insert only elements that it owns
 60:         locally (but any non-local elements will be sent to the
 61:         appropriate processor during matrix assembly).
 62:       - Always specify global rows and columns of matrix entries.
 63:    */
 64:   co1  = gamma * h * h / 2.0;
 65:   co2  = beta * h * h;
 66:   for (Ii=Istart; Ii<Iend; Ii++) {
 67:     i = Ii/(n2*n3); j = (Ii - i*n2*n3)/n3; k = Ii - i*n2*n3 - j*n3;
 68:     if (i>0) {
 69:       J    = Ii - n2*n3;  v = -1.0 + co1*(PetscScalar)i;
 70:       MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);
 71:     }
 72:     if (i<n1-1) {
 73:       J    = Ii + n2*n3;  v = -1.0 + co1*(PetscScalar)i;
 74:       MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);
 75:     }
 76:     if (j>0) {
 77:       J    = Ii - n3;  v = -1.0 + co1*(PetscScalar)j;
 78:       MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);
 79:     }
 80:     if (j<n2-1) {
 81:       J    = Ii + n3;  v = -1.0 + co1*(PetscScalar)j;
 82:       MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);
 83:     }
 84:     if (k>0) {
 85:       J    = Ii - 1;  v = -1.0 + co1*(PetscScalar)k;
 86:       MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);
 87:     }
 88:     if (k<n3-1) {
 89:       J    = Ii + 1;  v = -1.0 + co1*(PetscScalar)k;
 90:       MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);
 91:     }
 92:     v    = 6.0 + co2;
 93:     MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 94:   }
 95:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 96:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 98:   /* Create parallel vectors and Set right-hand side. */
 99:   VecCreate(PETSC_COMM_WORLD,&b);
100:   VecSetSizes(b,PETSC_DECIDE,n1*n2*n3);
101:   VecSetFromOptions(b);
102:   VecDuplicate(b,&x);
103:   VecDuplicate(b,&u);
104:   VecSet(b,1.0);

106:   /* Create linear solver context */
107:   KSPCreate(PETSC_COMM_WORLD,&ksp);
108:   KSPSetOperators(ksp,A,A);
109:   KSPSetTolerances(ksp,1.e-6,1.e-50,PETSC_DEFAULT,200);
110:   KSPSetFromOptions(ksp);

112:   /* Solve the linear system */
113:   KSPSolve(ksp,b,x);

115:   /* Free work space.  */
116:   KSPDestroy(&ksp);
117:   VecDestroy(&u);  VecDestroy(&x);
118:   VecDestroy(&b);  MatDestroy(&A);
119:   PetscFinalize();
120:   return 0;
121: }