Actual source code: ex64.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  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.


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

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

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


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

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

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

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

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

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

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

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


184: /*TEST

186:    test:
187:       nsize: 4
188:       requires: superlu_dist parmetis
189:       args: -ksp_monitor -pc_gasm_overlap 1 -sub_pc_type lu -sub_pc_factor_mat_solver_type superlu_dist -pc_gasm_total_subdomains 2

191: TEST*/