LCOV - code coverage report
Current view: top level - nep/tutorials - ex42.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 99 99 100.0 %
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             :    This example implements one of the problems found at
      12             :        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
      13             :        The University of Manchester.
      14             :    The details of the collection can be found at:
      15             :        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
      16             :            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
      17             : 
      18             :    The loaded_string problem is a rational eigenvalue problem for the
      19             :    finite element model of a loaded vibrating string.
      20             : */
      21             : 
      22             : static char help[] = "Illustrates computation of left eigenvectors and resolvent.\n\n"
      23             :   "This is based on loaded_string from the NLEVP collection.\n"
      24             :   "The command line options are:\n"
      25             :   "  -n <n>, dimension of the matrices.\n"
      26             :   "  -kappa <kappa>, stiffness of elastic spring.\n"
      27             :   "  -mass <m>, mass of the attached load.\n\n";
      28             : 
      29             : #include <slepcnep.h>
      30             : 
      31             : #define NMAT 3
      32             : 
      33           1 : int main(int argc,char **argv)
      34             : {
      35           1 :   Mat            A[NMAT];         /* problem matrices */
      36           1 :   FN             f[NMAT];         /* functions to define the nonlinear operator */
      37           1 :   NEP            nep;             /* nonlinear eigensolver context */
      38           1 :   RG             rg;
      39           1 :   Vec            v,r,z,w;
      40           1 :   PetscInt       n=100,Istart,Iend,i,nconv;
      41           1 :   PetscReal      kappa=1.0,m=1.0,nrm,tol;
      42           1 :   PetscScalar    lambda,sigma,numer[2],denom[2],omega1,omega2;
      43             : 
      44           1 :   PetscFunctionBeginUser;
      45           1 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      46             : 
      47           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      48           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
      49           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL));
      50           1 :   sigma = kappa/m;
      51           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string, n=%" PetscInt_FMT " kappa=%g m=%g\n\n",n,(double)kappa,(double)m));
      52             : 
      53             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      54             :                        Build the problem matrices
      55             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      56             : 
      57             :   /* initialize matrices */
      58           4 :   for (i=0;i<NMAT;i++) {
      59           3 :     PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
      60           3 :     PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
      61           3 :     PetscCall(MatSetFromOptions(A[i]));
      62             :   }
      63           1 :   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
      64             : 
      65             :   /* A0 */
      66         101 :   for (i=Istart;i<Iend;i++) {
      67         100 :     PetscCall(MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES));
      68         100 :     if (i>0) PetscCall(MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES));
      69         100 :     if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES));
      70             :   }
      71             : 
      72             :   /* A1 */
      73         101 :   for (i=Istart;i<Iend;i++) {
      74         100 :     PetscCall(MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES));
      75         100 :     if (i>0) PetscCall(MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES));
      76         100 :     if (i<n-1) PetscCall(MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES));
      77             :   }
      78             : 
      79             :   /* A2 */
      80           1 :   if (Istart<=n-1 && n-1<Iend) PetscCall(MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES));
      81             : 
      82             :   /* assemble matrices */
      83           4 :   for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
      84           4 :   for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
      85             : 
      86             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      87             :                        Create the problem functions
      88             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      89             : 
      90             :   /* f1=1 */
      91           1 :   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
      92           1 :   PetscCall(FNSetType(f[0],FNRATIONAL));
      93           1 :   numer[0] = 1.0;
      94           1 :   PetscCall(FNRationalSetNumerator(f[0],1,numer));
      95             : 
      96             :   /* f2=-lambda */
      97           1 :   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
      98           1 :   PetscCall(FNSetType(f[1],FNRATIONAL));
      99           1 :   numer[0] = -1.0; numer[1] = 0.0;
     100           1 :   PetscCall(FNRationalSetNumerator(f[1],2,numer));
     101             : 
     102             :   /* f3=lambda/(lambda-sigma) */
     103           1 :   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
     104           1 :   PetscCall(FNSetType(f[2],FNRATIONAL));
     105           1 :   numer[0] = 1.0; numer[1] = 0.0;
     106           1 :   denom[0] = 1.0; denom[1] = -sigma;
     107           1 :   PetscCall(FNRationalSetNumerator(f[2],2,numer));
     108           1 :   PetscCall(FNRationalSetDenominator(f[2],2,denom));
     109             : 
     110             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     111             :                 Create the eigensolver and solve the problem
     112             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     113             : 
     114           1 :   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
     115           1 :   PetscCall(NEPSetSplitOperator(nep,3,A,f,SUBSET_NONZERO_PATTERN));
     116           1 :   PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));
     117           1 :   PetscCall(NEPSetDimensions(nep,8,PETSC_DETERMINE,PETSC_DETERMINE));
     118             : 
     119             :   /* set two-sided NLEIGS solver */
     120           1 :   PetscCall(NEPSetType(nep,NEPNLEIGS));
     121           1 :   PetscCall(NEPNLEIGSSetFullBasis(nep,PETSC_TRUE));
     122           1 :   PetscCall(NEPSetTwoSided(nep,PETSC_TRUE));
     123           1 :   PetscCall(NEPGetRG(nep,&rg));
     124           1 :   PetscCall(RGSetType(rg,RGINTERVAL));
     125             : #if defined(PETSC_USE_COMPLEX)
     126             :   PetscCall(RGIntervalSetEndpoints(rg,4.0,700.0,-0.001,0.001));
     127             : #else
     128           1 :   PetscCall(RGIntervalSetEndpoints(rg,4.0,700.0,0,0));
     129             : #endif
     130           1 :   PetscCall(NEPSetTarget(nep,5.0));
     131             : 
     132           1 :   PetscCall(NEPSetFromOptions(nep));
     133           1 :   PetscCall(NEPSolve(nep));
     134             : 
     135             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     136             :                        Check left residual
     137             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     138           1 :   PetscCall(MatCreateVecs(A[0],&v,&r));
     139           1 :   PetscCall(VecDuplicate(v,&w));
     140           1 :   PetscCall(VecDuplicate(v,&z));
     141           1 :   PetscCall(NEPGetConverged(nep,&nconv));
     142           1 :   PetscCall(NEPGetTolerances(nep,&tol,NULL));
     143           9 :   for (i=0;i<nconv;i++) {
     144           8 :     PetscCall(NEPGetEigenpair(nep,i,&lambda,NULL,NULL,NULL));
     145           8 :     PetscCall(NEPGetLeftEigenvector(nep,i,v,NULL));
     146           8 :     PetscCall(NEPApplyAdjoint(nep,lambda,v,w,r,NULL,NULL));
     147           8 :     PetscCall(VecNorm(r,NORM_2,&nrm));
     148           8 :     if (nrm>tol*PetscAbsScalar(lambda)) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Left residual i=%" PetscInt_FMT " is above tolerance --> %g\n",i,(double)(nrm/PetscAbsScalar(lambda))));
     149             :   }
     150             : 
     151             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     152             :                       Operate with resolvent
     153             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     154           1 :   omega1 = 20.0;
     155           1 :   omega2 = 150.0;
     156           1 :   PetscCall(VecSet(v,0.0));
     157           1 :   PetscCall(VecSetValue(v,0,-1.0,INSERT_VALUES));
     158           1 :   PetscCall(VecSetValue(v,1,3.0,INSERT_VALUES));
     159           1 :   PetscCall(VecAssemblyBegin(v));
     160           1 :   PetscCall(VecAssemblyEnd(v));
     161           1 :   PetscCall(NEPApplyResolvent(nep,NULL,omega1,v,r));
     162           1 :   PetscCall(VecNorm(r,NORM_2,&nrm));
     163           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega1),(double)nrm));
     164           1 :   PetscCall(NEPApplyResolvent(nep,NULL,omega2,v,r));
     165           1 :   PetscCall(VecNorm(r,NORM_2,&nrm));
     166           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega2),(double)nrm));
     167           1 :   PetscCall(VecSet(v,1.0));
     168           1 :   PetscCall(NEPApplyResolvent(nep,NULL,omega1,v,r));
     169           1 :   PetscCall(VecNorm(r,NORM_2,&nrm));
     170           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega1),(double)nrm));
     171           1 :   PetscCall(NEPApplyResolvent(nep,NULL,omega2,v,r));
     172           1 :   PetscCall(VecNorm(r,NORM_2,&nrm));
     173           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega2),(double)nrm));
     174             : 
     175             :   /* clean up */
     176           1 :   PetscCall(NEPDestroy(&nep));
     177           4 :   for (i=0;i<NMAT;i++) {
     178           3 :     PetscCall(MatDestroy(&A[i]));
     179           3 :     PetscCall(FNDestroy(&f[i]));
     180             :   }
     181           1 :   PetscCall(VecDestroy(&v));
     182           1 :   PetscCall(VecDestroy(&r));
     183           1 :   PetscCall(VecDestroy(&w));
     184           1 :   PetscCall(VecDestroy(&z));
     185           1 :   PetscCall(SlepcFinalize());
     186             :   return 0;
     187             : }
     188             : 
     189             : /*TEST
     190             : 
     191             :    test:
     192             :       suffix: 1
     193             :       requires: !single
     194             : 
     195             : TEST*/

Generated by: LCOV version 1.14