Actual source code: test5.c
slepc-3.17.0 2022-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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[] = "Test matrix rational function.\n\n";
13: #include <slepcfn.h>
15: /*
16: Compute matrix rational function B = q(A)\p(A)
17: */
18: PetscErrorCode TestMatRational(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
19: {
20: PetscBool set,flg;
21: PetscInt n;
22: Mat F;
23: Vec v,f0;
24: PetscReal nrm;
27: MatGetSize(A,&n,NULL);
28: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
29: PetscObjectSetName((PetscObject)F,"F");
30: /* compute matrix function */
31: if (inplace) {
32: MatCopy(A,F,SAME_NONZERO_PATTERN);
33: MatIsHermitianKnown(A,&set,&flg);
34: if (set && flg) MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE);
35: FNEvaluateFunctionMat(fn,F,NULL);
36: } else FNEvaluateFunctionMat(fn,A,F);
37: if (verbose) {
38: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
39: MatView(A,viewer);
40: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
41: MatView(F,viewer);
42: }
43: /* print matrix norm for checking */
44: MatNorm(F,NORM_1,&nrm);
45: PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm);
46: /* check FNEvaluateFunctionMatVec() */
47: MatCreateVecs(A,&v,&f0);
48: MatGetColumnVector(F,f0,0);
49: FNEvaluateFunctionMatVec(fn,A,v);
50: VecAXPY(v,-1.0,f0);
51: VecNorm(v,NORM_2,&nrm);
52: if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
53: MatDestroy(&F);
54: VecDestroy(&v);
55: VecDestroy(&f0);
56: PetscFunctionReturn(0);
57: }
59: int main(int argc,char **argv)
60: {
61: FN fn;
62: Mat A;
63: PetscInt i,j,n=10,np,nq;
64: PetscScalar *As,p[10],q[10];
65: PetscViewer viewer;
66: PetscBool verbose,inplace;
68: SlepcInitialize(&argc,&argv,(char*)0,help);
69: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
70: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
71: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
72: PetscPrintf(PETSC_COMM_WORLD,"Matrix rational function, n=%" PetscInt_FMT ".\n",n);
74: /* Create rational function r(x)=p(x)/q(x) */
75: FNCreate(PETSC_COMM_WORLD,&fn);
76: FNSetType(fn,FNRATIONAL);
77: np = 2; nq = 3;
78: p[0] = -3.1; p[1] = 1.1;
79: q[0] = 1.0; q[1] = -2.0; q[2] = 3.5;
80: FNRationalSetNumerator(fn,np,p);
81: FNRationalSetDenominator(fn,nq,q);
82: FNSetFromOptions(fn);
84: /* Set up viewer */
85: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
86: FNView(fn,viewer);
87: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
89: /* Create matrices */
90: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
91: PetscObjectSetName((PetscObject)A,"A");
93: /* Fill A with a symmetric Toeplitz matrix */
94: MatDenseGetArray(A,&As);
95: for (i=0;i<n;i++) As[i+i*n]=2.0;
96: for (j=1;j<3;j++) {
97: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
98: }
99: MatDenseRestoreArray(A,&As);
100: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
101: TestMatRational(fn,A,viewer,verbose,inplace);
103: /* Repeat with same matrix as non-symmetric */
104: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
105: TestMatRational(fn,A,viewer,verbose,inplace);
107: MatDestroy(&A);
108: FNDestroy(&fn);
109: SlepcFinalize();
110: return 0;
111: }
113: /*TEST
115: testset:
116: output_file: output/test5_1.out
117: requires: !single
118: test:
119: suffix: 1
120: test:
121: suffix: 2
122: args: -inplace
124: TEST*/