Actual source code: ex64.c

  1: /*
  2:  *
  3:  *  Created on: Aug 10, 2015
  4:  *      Author: Fande Kong  <fdkong.jd@gmail.com>
  5:  */

  7: static char help[] = "Illustrates use of the preconditioner GASM.\n \
  8:    using hierarchical partitioning and MatIncreaseOverlapSplit \
  9:         -pc_gasm_total_subdomains\n \
 10:         -pc_gasm_print_subdomains\n \n";

 12: /*
 13:    Note:  This example focuses on setting the subdomains for the GASM
 14:    preconditioner for a problem on a 2D rectangular grid.  See ex1.c
 15:    and ex2.c for more detailed comments on the basic usage of KSP
 16:    (including working with matrices and vectors).

 18:    The GASM preconditioner is fully parallel.  The user-space routine
 19:    CreateSubdomains2D that computes the domain decomposition is also parallel
 20:    and attempts to generate both subdomains straddling processors and multiple
 21:    domains per processor.

 23:    This matrix in this linear system arises from the discretized Laplacian,
 24:    and thus is not very interesting in terms of experimenting with variants
 25:    of the GASM preconditioner.
 26: */

 28: /*T
 29:    Concepts: KSP^Additive Schwarz Method (GASM) with user-defined subdomains
 30:    Processors: n
 31: T*/

 33: /*
 34:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 35:   automatically includes:
 36:      petscsys.h    - base PETSc routines   petscvec.h - vectors
 37:      petscmat.h    - matrices
 38:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 39:      petscviewer.h - viewers               petscpc.h  - preconditioners
 40: */
 41: #include <petscksp.h>
 42: #include <petscmat.h>

 44: int main(int argc,char **args)
 45: {
 46:   Vec            x,b,u;                  /* approx solution, RHS, exact solution */
 47:   Mat            A,perA;                      /* linear system matrix */
 48:   KSP            ksp;                    /* linear solver context */
 49:   PC             pc;                     /* PC context */
 50:   PetscInt       overlap;                /* width of subdomain overlap */
 51:   PetscInt       m,n;                    /* mesh dimensions in x- and y- directions */
 52:   PetscInt       M,N;                    /* number of subdomains in x- and y- directions */
 53:   PetscInt       i,j,Ii,J,Istart,Iend;
 55:   PetscMPIInt    size;
 56:   PetscBool      flg;
 57:   PetscBool       user_set_subdomains;
 58:   PetscScalar     v, one = 1.0;
 59:   MatPartitioning part;
 60:   IS              coarseparts,fineparts;
 61:   IS              is,isn,isrows;
 62:   MPI_Comm        comm;
 63:   PetscReal       e;

 65:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 66:   comm = PETSC_COMM_WORLD;
 67:   MPI_Comm_size(comm,&size);
 68:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex62","PC");
 69:   m = 15;
 70:   PetscOptionsInt("-M", "Number of mesh points in the x-direction","PCGASMCreateSubdomains2D",m,&m,NULL);
 71:   n = 17;
 72:   PetscOptionsInt("-N","Number of mesh points in the y-direction","PCGASMCreateSubdomains2D",n,&n,NULL);
 73:   user_set_subdomains = PETSC_FALSE;
 74:   PetscOptionsBool("-user_set_subdomains","Use the user-specified 2D tiling of mesh by subdomains","PCGASMCreateSubdomains2D",user_set_subdomains,&user_set_subdomains,NULL);
 75:   M = 2;
 76:   PetscOptionsInt("-Mdomains","Number of subdomain tiles in the x-direction","PCGASMSetSubdomains2D",M,&M,NULL);
 77:   N = 1;
 78:   PetscOptionsInt("-Ndomains","Number of subdomain tiles in the y-direction","PCGASMSetSubdomains2D",N,&N,NULL);
 79:   overlap = 1;
 80:   PetscOptionsInt("-overlap","Size of tile overlap.","PCGASMSetSubdomains2D",overlap,&overlap,NULL);
 81:   PetscOptionsEnd();

 83:   /* -------------------------------------------------------------------
 84:          Compute the matrix and right-hand-side vector that define
 85:          the linear system, Ax = b.
 86:      ------------------------------------------------------------------- */

 88:   /*
 89:      Assemble the matrix for the five point stencil, YET AGAIN
 90:   */
 91:   MatCreate(comm,&A);
 92:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 93:   MatSetFromOptions(A);
 94:   MatSetUp(A);
 95:   MatGetOwnershipRange(A,&Istart,&Iend);
 96:   for (Ii=Istart; Ii<Iend; Ii++) {
 97:     v = -1.0; i = Ii/n; j = Ii - i*n;
 98:     if (i>0)   {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 99:     if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
100:     if (j>0)   {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
101:     if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
102:     v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
103:   }
104:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

107:   /*
108:     Partition the graph of the matrix
109:   */
110:   MatPartitioningCreate(comm,&part);
111:   MatPartitioningSetAdjacency(part,A);
112:   MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);
113:   MatPartitioningHierarchicalSetNcoarseparts(part,2);
114:   MatPartitioningHierarchicalSetNfineparts(part,2);
115:   MatPartitioningSetFromOptions(part);
116:   /* get new processor owner number of each vertex */
117:   MatPartitioningApply(part,&is);
118:   /* get coarse parts according to which we rearrange the matrix */
119:   MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts);
120:   /* fine parts are used with coarse parts */
121:   MatPartitioningHierarchicalGetFineparts(part,&fineparts);
122:   /* get new global number of each old global number */
123:   ISPartitioningToNumbering(is,&isn);
124:   ISBuildTwoSided(is,NULL,&isrows);
125:   MatCreateSubMatrix(A,isrows,isrows,MAT_INITIAL_MATRIX,&perA);
126:   MatCreateVecs(perA,&b,NULL);
127:   VecSetFromOptions(b);
128:   VecDuplicate(b,&u);
129:   VecDuplicate(b,&x);
130:   VecSet(u,one);
131:   MatMult(perA,u,b);
132:   KSPCreate(comm,&ksp);
133:   /*
134:      Set operators. Here the matrix that defines the linear system
135:      also serves as the preconditioning matrix.
136:   */
137:   KSPSetOperators(ksp,perA,perA);

139:   /*
140:      Set the default preconditioner for this program to be GASM
141:   */
142:   KSPGetPC(ksp,&pc);
143:   PCSetType(pc,PCGASM);
144:   KSPSetFromOptions(ksp);
145:   /* -------------------------------------------------------------------
146:                       Solve the linear system
147:      ------------------------------------------------------------------- */

149:   KSPSolve(ksp,b,x);
150:   /* -------------------------------------------------------------------
151:                       Compare result to the exact solution
152:      ------------------------------------------------------------------- */
153:   VecAXPY(x,-1.0,u);
154:   VecNorm(x,NORM_INFINITY, &e);

156:   flg  = PETSC_FALSE;
157:   PetscOptionsGetBool(NULL,NULL,"-print_error",&flg,NULL);
158:   if (flg) {
159:     PetscPrintf(PETSC_COMM_WORLD, "Infinity norm of the error: %g\n", (double)e);
160:   }
161:   /*
162:      Free work space.  All PETSc objects should be destroyed when they
163:      are no longer needed.
164:   */
165:   KSPDestroy(&ksp);
166:   VecDestroy(&u);
167:   VecDestroy(&x);
168:   VecDestroy(&b);
169:   MatDestroy(&A);
170:   MatDestroy(&perA);
171:   ISDestroy(&is);
172:   ISDestroy(&coarseparts);
173:   ISDestroy(&fineparts);
174:   ISDestroy(&isrows);
175:   ISDestroy(&isn);
176:   MatPartitioningDestroy(&part);
177:   PetscFinalize();
178:   return ierr;
179: }

181: /*TEST

183:    test:
184:       nsize: 4
185:       requires: superlu_dist parmetis
186:       args: -ksp_monitor -pc_gasm_overlap 1 -sub_pc_type lu -sub_pc_factor_mat_solver_type superlu_dist -pc_gasm_total_subdomains 2

188: TEST*/