Actual source code: test11.c
slepc-3.17.1 2022-04-11
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: */
10: /*
11: Define the function
13: f(x) = (exp(x)-1)/x (the phi_1 function)
15: with the following tree:
17: f(x) f(x) (combined by division)
18: / \ p(x) = x (polynomial)
19: a(x) p(x) a(x) (combined by addition)
20: / \ e(x) = exp(x) (exponential)
21: e(x) c(x) c(x) = -1 (constant)
22: */
24: static char help[] = "Another test of a combined function.\n\n";
26: #include <slepcfn.h>
28: /*
29: Compute matrix function B = A\(exp(A)-I)
30: */
31: PetscErrorCode TestMatCombine(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
32: {
33: PetscBool set,flg;
34: PetscInt n;
35: Mat F;
36: Vec v,f0;
37: PetscReal nrm;
40: MatGetSize(A,&n,NULL);
41: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
42: PetscObjectSetName((PetscObject)F,"F");
43: /* compute matrix function */
44: if (inplace) {
45: MatCopy(A,F,SAME_NONZERO_PATTERN);
46: MatIsHermitianKnown(A,&set,&flg);
47: if (set && flg) MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE);
48: FNEvaluateFunctionMat(fn,F,NULL);
49: } else FNEvaluateFunctionMat(fn,A,F);
50: if (verbose) {
51: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
52: MatView(A,viewer);
53: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
54: MatView(F,viewer);
55: }
56: /* print matrix norm for checking */
57: MatNorm(F,NORM_1,&nrm);
58: PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %6.3f\n",(double)nrm);
59: /* check FNEvaluateFunctionMatVec() */
60: MatCreateVecs(A,&v,&f0);
61: MatGetColumnVector(F,f0,0);
62: FNEvaluateFunctionMatVec(fn,A,v);
63: VecAXPY(v,-1.0,f0);
64: VecNorm(v,NORM_2,&nrm);
65: if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
66: MatDestroy(&F);
67: VecDestroy(&v);
68: VecDestroy(&f0);
69: PetscFunctionReturn(0);
70: }
72: int main(int argc,char **argv)
73: {
74: FN f,p,a,e,c,f1,f2;
75: FNCombineType ctype;
76: Mat A;
77: PetscInt i,j,n=10,np;
78: PetscScalar x,y,yp,*As,coeffs[10];
79: char strx[50],str[50];
80: PetscViewer viewer;
81: PetscBool verbose,inplace;
83: SlepcInitialize(&argc,&argv,(char*)0,help);
84: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
85: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
86: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
87: PetscPrintf(PETSC_COMM_WORLD,"Phi1 via a combined function, n=%" PetscInt_FMT ".\n",n);
89: /* Create function */
91: /* e(x) = exp(x) */
92: FNCreate(PETSC_COMM_WORLD,&e);
93: PetscObjectSetName((PetscObject)e,"e");
94: FNSetType(e,FNEXP);
95: FNSetFromOptions(e);
96: /* c(x) = -1 */
97: FNCreate(PETSC_COMM_WORLD,&c);
98: PetscObjectSetName((PetscObject)c,"c");
99: FNSetType(c,FNRATIONAL);
100: FNSetFromOptions(c);
101: np = 1;
102: coeffs[0] = -1.0;
103: FNRationalSetNumerator(c,np,coeffs);
104: /* a(x) */
105: FNCreate(PETSC_COMM_WORLD,&a);
106: PetscObjectSetName((PetscObject)a,"a");
107: FNSetType(a,FNCOMBINE);
108: FNSetFromOptions(a);
109: FNCombineSetChildren(a,FN_COMBINE_ADD,e,c);
110: /* p(x) = x */
111: FNCreate(PETSC_COMM_WORLD,&p);
112: PetscObjectSetName((PetscObject)p,"p");
113: FNSetType(p,FNRATIONAL);
114: FNSetFromOptions(p);
115: np = 2;
116: coeffs[0] = 1.0; coeffs[1] = 0.0;
117: FNRationalSetNumerator(p,np,coeffs);
118: /* f(x) */
119: FNCreate(PETSC_COMM_WORLD,&f);
120: PetscObjectSetName((PetscObject)f,"f");
121: FNSetType(f,FNCOMBINE);
122: FNSetFromOptions(f);
123: FNCombineSetChildren(f,FN_COMBINE_DIVIDE,a,p);
125: /* Set up viewer */
126: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
127: FNCombineGetChildren(f,&ctype,&f1,&f2);
128: PetscPrintf(PETSC_COMM_WORLD,"Two functions combined with division:\n");
129: FNView(f1,viewer);
130: FNView(f2,viewer);
131: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
133: /* Scalar evaluation */
134: x = 2.2;
135: SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
136: FNEvaluateFunction(f,x,&y);
137: FNEvaluateDerivative(f,x,&yp);
138: SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
139: PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str);
140: SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
141: PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str);
143: /* Create matrices */
144: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
145: PetscObjectSetName((PetscObject)A,"A");
147: /* Fill A with 1-D Laplacian matrix */
148: MatDenseGetArray(A,&As);
149: for (i=0;i<n;i++) As[i+i*n]=2.0;
150: j=1;
151: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=-1.0; As[(i+j)+i*n]=-1.0; }
152: MatDenseRestoreArray(A,&As);
153: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
154: TestMatCombine(f,A,viewer,verbose,inplace);
156: /* Repeat with same matrix as non-symmetric */
157: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
158: TestMatCombine(f,A,viewer,verbose,inplace);
160: MatDestroy(&A);
161: FNDestroy(&f);
162: FNDestroy(&p);
163: FNDestroy(&a);
164: FNDestroy(&e);
165: FNDestroy(&c);
166: SlepcFinalize();
167: return 0;
168: }
170: /*TEST
172: test:
173: suffix: 1
174: nsize: 1
176: test:
177: suffix: 2
178: nsize: 1
179: args: -inplace
180: output_file: output/test11_1.out
182: TEST*/