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