Actual source code: test6.c
slepc-3.7.3 2016-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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: */
21: /*
22: Define the function
24: f(x) = (1-x^2) exp( -x/(1+x^2) )
26: with the following tree:
28: f(x) f(x) (combined by product)
29: / \ g(x) = 1-x^2 (polynomial)
30: g(x) h(x) h(x) (combined by composition)
31: / \ r(x) = -x/(1+x^2) (rational)
32: r(x) e(x) e(x) = exp(x) (exponential)
33: */
35: static char help[] = "Test combined function.\n\n";
37: #include <slepcfn.h>
41: /*
42: Compute matrix function B = (I-A^2) exp( -(I+A^2)\A )
43: */
44: PetscErrorCode TestMatCombine(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
45: {
47: PetscBool set,flg;
48: PetscInt n;
49: Mat F;
50: Vec v,f0;
51: PetscReal nrm;
54: MatGetSize(A,&n,NULL);
55: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
56: PetscObjectSetName((PetscObject)F,"F");
57: /* compute square root */
58: if (inplace) {
59: MatCopy(A,F,SAME_NONZERO_PATTERN);
60: MatIsHermitianKnown(A,&set,&flg);
61: if (set && flg) { MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE); }
62: FNEvaluateFunctionMat(fn,F,NULL);
63: } else {
64: FNEvaluateFunctionMat(fn,A,F);
65: }
66: if (verbose) {
67: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
68: MatView(A,viewer);
69: PetscPrintf(PETSC_COMM_WORLD,"Computed expm(A) - - - - - - -\n");
70: MatView(F,viewer);
71: }
72: /* print matrix norm for checking */
73: MatNorm(F,NORM_1,&nrm);
74: PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm);
75: /* check FNEvaluateFunctionMatVec() */
76: MatCreateVecs(A,&v,&f0);
77: MatGetColumnVector(F,f0,0);
78: FNEvaluateFunctionMatVec(fn,A,v);
79: VecAXPY(v,-1.0,f0);
80: VecNorm(v,NORM_2,&nrm);
81: if (nrm>100*PETSC_MACHINE_EPSILON) {
82: PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
83: }
84: MatDestroy(&F);
85: VecDestroy(&v);
86: VecDestroy(&f0);
87: return(0);
88: }
92: int main(int argc,char **argv)
93: {
95: FN f,g,h,e,r;
96: Mat A;
97: PetscInt i,j,n=10,np,nq;
98: PetscScalar x,y,yp,*As,p[10],q[10];
99: char strx[50],str[50];
100: PetscViewer viewer;
101: PetscBool verbose,inplace;
103: SlepcInitialize(&argc,&argv,(char*)0,help);
104: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
105: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
106: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
107: PetscPrintf(PETSC_COMM_WORLD,"Combined function, n=%D.\n",n);
109: /* Create function */
111: /* e(x) = exp(x) */
112: FNCreate(PETSC_COMM_WORLD,&e);
113: FNSetType(e,FNEXP);
114: /* r(x) = x/(1+x^2) */
115: FNCreate(PETSC_COMM_WORLD,&r);
116: FNSetType(r,FNRATIONAL);
117: np = 2; nq = 3;
118: p[0] = -1.0; p[1] = 0.0;
119: q[0] = 1.0; q[1] = 0.0; q[2] = 1.0;
120: FNRationalSetNumerator(r,np,p);
121: FNRationalSetDenominator(r,nq,q);
122: /* h(x) */
123: FNCreate(PETSC_COMM_WORLD,&h);
124: FNSetType(h,FNCOMBINE);
125: FNCombineSetChildren(h,FN_COMBINE_COMPOSE,r,e);
126: /* g(x) = 1-x^2 */
127: FNCreate(PETSC_COMM_WORLD,&g);
128: FNSetType(g,FNRATIONAL);
129: np = 3;
130: p[0] = -1.0; p[1] = 0.0; p[2] = 1.0;
131: FNRationalSetNumerator(g,np,p);
132: /* f(x) */
133: FNCreate(PETSC_COMM_WORLD,&f);
134: FNSetType(f,FNCOMBINE);
135: FNCombineSetChildren(f,FN_COMBINE_MULTIPLY,g,h);
137: /* Set up viewer */
138: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
139: FNView(f,viewer);
140: if (verbose) {
141: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
142: }
144: /* Scalar evaluation */
145: x = 2.2;
146: SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
147: FNEvaluateFunction(f,x,&y);
148: FNEvaluateDerivative(f,x,&yp);
149: SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
150: PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str);
151: SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
152: PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str);
154: /* Create matrices */
155: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
156: PetscObjectSetName((PetscObject)A,"A");
158: /* Fill A with a symmetric Toeplitz matrix */
159: MatDenseGetArray(A,&As);
160: for (i=0;i<n;i++) As[i+i*n]=2.0;
161: for (j=1;j<3;j++) {
162: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
163: }
164: MatDenseRestoreArray(A,&As);
165: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
166: TestMatCombine(f,A,viewer,verbose,inplace);
168: /* Repeat with same matrix as non-symmetric */
169: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
170: TestMatCombine(f,A,viewer,verbose,inplace);
172: MatDestroy(&A);
173: FNDestroy(&f);
174: FNDestroy(&g);
175: FNDestroy(&h);
176: FNDestroy(&e);
177: FNDestroy(&r);
178: SlepcFinalize();
179: return ierr;
180: }