LCOV - code coverage report
Current view: top level - mfn/tutorials - ex26.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 56 57 98.2 %
Date: 2024-12-18 00:42:09 Functions: 1 1 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : 
      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";
      16             : 
      17             : #include <slepcmfn.h>
      18             : 
      19           1 : int main(int argc,char **argv)
      20             : {
      21           1 :   Mat            A;           /* problem matrix */
      22           1 :   MFN            mfn;
      23           1 :   FN             f;
      24           1 :   PetscReal      norm,tol;
      25           1 :   Vec            v,y,z;
      26           1 :   PetscInt       N,n=10,m,Istart,Iend,i,j,II;
      27           1 :   PetscBool      flag;
      28             : 
      29           1 :   PetscFunctionBeginUser;
      30           1 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      31             : 
      32           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      33           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
      34           1 :   if (!flag) m=n;
      35           1 :   N = n*m;
      36           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
      37             : 
      38             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      39             :                  Compute the discrete 2-D Laplacian, A
      40             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      41             : 
      42           1 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
      43           1 :   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
      44           1 :   PetscCall(MatSetFromOptions(A));
      45             : 
      46           1 :   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
      47         101 :   for (II=Istart;II<Iend;II++) {
      48         100 :     i = II/n; j = II-i*n;
      49         100 :     if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
      50         100 :     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
      51         100 :     if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
      52         100 :     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
      53         100 :     PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
      54             :   }
      55             : 
      56           1 :   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
      57           1 :   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
      58             : 
      59             :   /* set symmetry flag so that solver can exploit it */
      60           1 :   PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
      61             : 
      62             :   /* set v = e_1 */
      63           1 :   PetscCall(MatCreateVecs(A,NULL,&v));
      64           1 :   PetscCall(VecSetValue(v,0,1.0,INSERT_VALUES));
      65           1 :   PetscCall(VecAssemblyBegin(v));
      66           1 :   PetscCall(VecAssemblyEnd(v));
      67           1 :   PetscCall(VecDuplicate(v,&y));
      68           1 :   PetscCall(VecDuplicate(v,&z));
      69             : 
      70             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      71             :              Create the solver, set the matrix and the function
      72             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      73           1 :   PetscCall(MFNCreate(PETSC_COMM_WORLD,&mfn));
      74           1 :   PetscCall(MFNSetOperator(mfn,A));
      75           1 :   PetscCall(MFNGetFN(mfn,&f));
      76           1 :   PetscCall(FNSetType(f,FNSQRT));
      77           1 :   PetscCall(MFNSetErrorIfNotConverged(mfn,PETSC_TRUE));
      78           1 :   PetscCall(MFNSetFromOptions(mfn));
      79             : 
      80             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      81             :                       First solve: y=sqrt(A)*v
      82             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      83             : 
      84           1 :   PetscCall(MFNSolve(mfn,v,y));
      85           1 :   PetscCall(VecNorm(y,NORM_2,&norm));
      86           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Intermediate vector has norm %g\n",(double)norm));
      87             : 
      88             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      89             :              Second solve: z=sqrt(A)*y and compare against A*v
      90             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      91             : 
      92           1 :   PetscCall(MFNSolve(mfn,y,z));
      93           1 :   PetscCall(MFNGetTolerances(mfn,&tol,NULL));
      94             : 
      95           1 :   PetscCall(MatMult(A,v,y));   /* overwrite y */
      96           1 :   PetscCall(VecAXPY(y,-1.0,z));
      97           1 :   PetscCall(VecNorm(y,NORM_2,&norm));
      98             : 
      99           1 :   if (norm<tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error norm is less than the requested tolerance\n\n"));
     100           0 :   else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error norm larger than tolerance: %3.1e\n\n",(double)norm));
     101             : 
     102             :   /*
     103             :      Free work space
     104             :   */
     105           1 :   PetscCall(MFNDestroy(&mfn));
     106           1 :   PetscCall(MatDestroy(&A));
     107           1 :   PetscCall(VecDestroy(&v));
     108           1 :   PetscCall(VecDestroy(&y));
     109           1 :   PetscCall(VecDestroy(&z));
     110           1 :   PetscCall(SlepcFinalize());
     111             :   return 0;
     112             : }
     113             : 
     114             : /*TEST
     115             : 
     116             :    test:
     117             :       suffix: 1
     118             :       args: -mfn_tol 1e-4
     119             : 
     120             : TEST*/

Generated by: LCOV version 1.14