Actual source code: ex26.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Computes the action of the square root of the 2-D Laplacian.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n"
 15:   "To draw the solution run with -mfn_view_solution draw -draw_pause -1\n\n";

 17: #include <slepcmfn.h>

 19: int main(int argc,char **argv)
 20: {
 21:   Mat            A;           /* problem matrix */
 22:   MFN            mfn;
 23:   FN             f;
 24:   PetscReal      norm,tol;
 25:   Vec            v,y,z;
 26:   PetscInt       N,n=10,m,Istart,Iend,i,j,II;
 28:   PetscBool      flag;

 30:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 32:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 33:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 34:   if (!flag) m=n;
 35:   N = n*m;
 36:   PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%D (%Dx%D grid)\n\n",N,n,m);

 38:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 39:                  Compute the discrete 2-D Laplacian, A
 40:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 42:   MatCreate(PETSC_COMM_WORLD,&A);
 43:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 44:   MatSetFromOptions(A);
 45:   MatSetUp(A);

 47:   MatGetOwnershipRange(A,&Istart,&Iend);
 48:   for (II=Istart;II<Iend;II++) {
 49:     i = II/n; j = II-i*n;
 50:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 51:     if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 52:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 53:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 54:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 55:   }

 57:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 58:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 60:   /* set symmetry flag so that solver can exploit it */
 61:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

 63:   /* set v = e_1 */
 64:   MatCreateVecs(A,NULL,&v);
 65:   VecSetValue(v,0,1.0,INSERT_VALUES);
 66:   VecAssemblyBegin(v);
 67:   VecAssemblyEnd(v);
 68:   VecDuplicate(v,&y);
 69:   VecDuplicate(v,&z);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72:              Create the solver, set the matrix and the function
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 74:   MFNCreate(PETSC_COMM_WORLD,&mfn);
 75:   MFNSetOperator(mfn,A);
 76:   MFNGetFN(mfn,&f);
 77:   FNSetType(f,FNSQRT);
 78:   MFNSetErrorIfNotConverged(mfn,PETSC_TRUE);
 79:   MFNSetFromOptions(mfn);

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:                       First solve: y=sqrt(A)*v
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   MFNSolve(mfn,v,y);
 86:   VecNorm(y,NORM_2,&norm);
 87:   PetscPrintf(PETSC_COMM_WORLD," Intermediate vector has norm %g\n",(double)norm);

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 90:              Second solve: z=sqrt(A)*y and compare against A*v
 91:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 93:   MFNSolve(mfn,y,z);
 94:   MFNGetTolerances(mfn,&tol,NULL);

 96:   MatMult(A,v,y);   /* overwrite y */
 97:   VecAXPY(y,-1.0,z);
 98:   VecNorm(y,NORM_2,&norm);

100:   if (norm<tol) {
101:     PetscPrintf(PETSC_COMM_WORLD," Error norm is less than the requested tolerance\n\n");
102:   } else {
103:     PetscPrintf(PETSC_COMM_WORLD," Error norm larger than tolerance: %3.1e\n\n",(double)norm);
104:   }

106:   /*
107:      Free work space
108:   */
109:   MFNDestroy(&mfn);
110:   MatDestroy(&A);
111:   VecDestroy(&v);
112:   VecDestroy(&y);
113:   VecDestroy(&z);
114:   SlepcFinalize();
115:   return ierr;
116: }

118: /*TEST

120:    test:
121:       suffix: 1
122:       args: -mfn_tol 1e-4

124: TEST*/