Actual source code: test5.c

slepc-3.6.3 2016-03-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Test matrix rational function.\n\n";

 24: #include <slepcfn.h>

 28: int main(int argc,char **argv)
 29: {
 31:   FN             fn;
 32:   Mat            A,B;
 33:   PetscInt       i,j,n=10,np,nq;
 34:   PetscReal      nrm;
 35:   PetscScalar    *As,p[10],q[10];
 36:   PetscViewer    viewer;
 37:   PetscBool      verbose;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 41:   PetscOptionsHasName(NULL,"-verbose",&verbose);
 42:   PetscPrintf(PETSC_COMM_WORLD,"Matrix rational function, n=%D.\n",n);

 44:   /* Create rational function r(x)=p(x)/q(x) */
 45:   FNCreate(PETSC_COMM_WORLD,&fn);
 46:   FNSetType(fn,FNRATIONAL);
 47:   np = 2; nq = 3;
 48:   p[0] = -3.1; p[1] = 1.1;
 49:   q[0] = 1.0; q[1] = -2.0; q[2] = 3.5;
 50:   FNRationalSetNumerator(fn,np,p);
 51:   FNRationalSetDenominator(fn,nq,q);
 52:   FNSetFromOptions(fn);

 54:   /* Set up viewer */
 55:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 56:   FNView(fn,viewer);
 57:   if (verbose) {
 58:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 59:   }

 61:   /* Create matrices */
 62:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
 63:   PetscObjectSetName((PetscObject)A,"A");
 64:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&B);
 65:   PetscObjectSetName((PetscObject)B,"B");

 67:   /* Fill A with a symmetric Toeplitz matrix */
 68:   MatDenseGetArray(A,&As);
 69:   for (i=0;i<n;i++) As[i+i*n]=2.0;
 70:   for (j=1;j<3;j++) {
 71:     for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
 72:   }
 73:   MatDenseRestoreArray(A,&As);
 74:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
 75:   if (verbose) {
 76:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 77:     MatView(A,viewer);
 78:   }

 80:   /* Evaluate matrix function */
 81:   FNEvaluateFunctionMat(fn,A,B);
 82:   if (verbose) {
 83:     PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
 84:     MatView(B,viewer);
 85:   }
 86:   MatNorm(B,NORM_1,&nrm);
 87:   PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm);

 89:   /* Repeat with same matrix as non-symmetric */
 90:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);

 92:   /* Evaluate matrix function */
 93:   FNEvaluateFunctionMat(fn,A,B);
 94:   if (verbose) {
 95:     PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
 96:     MatView(B,viewer);
 97:   }
 98:   MatNorm(B,NORM_1,&nrm);
 99:   PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm);

101:   MatDestroy(&A);
102:   MatDestroy(&B);
103:   FNDestroy(&fn);
104:   SlepcFinalize();
105:   return 0;
106: }