Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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
12 :
13 : f(x) = (1-x^2) exp(-x/(1+x^2))
14 :
15 : with the following tree:
16 :
17 : f(x) f(x) (combined by product)
18 : / \ g(x) = 1-x^2 (polynomial)
19 : g(x) h(x) h(x) (combined by composition)
20 : / \ r(x) = -x/(1+x^2) (rational)
21 : r(x) e(x) e(x) = exp(x) (exponential)
22 : */
23 :
24 : static char help[] = "Test combined function.\n\n";
25 :
26 : #include <slepcfn.h>
27 :
28 : /*
29 : Compute matrix function B = (I-A^2) exp(-(I+A^2)\A)
30 : */
31 4 : PetscErrorCode TestMatCombine(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
32 : {
33 4 : PetscBool set,flg;
34 4 : PetscInt n;
35 4 : Mat F,Acopy;
36 4 : Vec v,f0;
37 4 : PetscReal nrm;
38 :
39 4 : PetscFunctionBeginUser;
40 4 : PetscCall(MatGetSize(A,&n,NULL));
41 4 : PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&F));
42 4 : PetscCall(PetscObjectSetName((PetscObject)F,"F"));
43 : /* compute matrix function */
44 4 : if (inplace) {
45 2 : PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN));
46 2 : PetscCall(MatIsHermitianKnown(A,&set,&flg));
47 2 : if (set && flg) PetscCall(MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE));
48 2 : PetscCall(FNEvaluateFunctionMat(fn,F,NULL));
49 : } else {
50 2 : PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Acopy));
51 2 : PetscCall(FNEvaluateFunctionMat(fn,A,F));
52 : /* check that A has not been modified */
53 2 : PetscCall(MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN));
54 2 : PetscCall(MatNorm(Acopy,NORM_1,&nrm));
55 2 : if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm));
56 2 : PetscCall(MatDestroy(&Acopy));
57 : }
58 4 : if (verbose) {
59 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n"));
60 0 : PetscCall(MatView(A,viewer));
61 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n"));
62 0 : PetscCall(MatView(F,viewer));
63 : }
64 : /* print matrix norm for checking */
65 4 : PetscCall(MatNorm(F,NORM_1,&nrm));
66 4 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %6.3f\n",(double)nrm));
67 : /* check FNEvaluateFunctionMatVec() */
68 4 : PetscCall(MatCreateVecs(A,&v,&f0));
69 4 : PetscCall(MatGetColumnVector(F,f0,0));
70 4 : PetscCall(FNEvaluateFunctionMatVec(fn,A,v));
71 4 : PetscCall(VecAXPY(v,-1.0,f0));
72 4 : PetscCall(VecNorm(v,NORM_2,&nrm));
73 4 : if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm));
74 4 : PetscCall(MatDestroy(&F));
75 4 : PetscCall(VecDestroy(&v));
76 4 : PetscCall(VecDestroy(&f0));
77 4 : PetscFunctionReturn(PETSC_SUCCESS);
78 : }
79 :
80 2 : int main(int argc,char **argv)
81 : {
82 2 : FN f,g,h,e,r,fcopy;
83 2 : Mat A=NULL;
84 2 : PetscInt i,j,n=10,np,nq;
85 2 : PetscScalar x,y,yp,*As,p[10],q[10];
86 2 : char strx[50],str[50];
87 2 : PetscViewer viewer;
88 2 : PetscBool verbose,inplace,matcuda;
89 :
90 2 : PetscFunctionBeginUser;
91 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
92 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
93 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
94 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-inplace",&inplace));
95 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda));
96 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Combined function, n=%" PetscInt_FMT ".\n",n));
97 :
98 : /* Create function */
99 :
100 : /* e(x) = exp(x) */
101 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&e));
102 2 : PetscCall(FNSetType(e,FNEXP));
103 2 : PetscCall(FNSetFromOptions(e));
104 : /* r(x) = x/(1+x^2) */
105 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&r));
106 2 : PetscCall(FNSetType(r,FNRATIONAL));
107 2 : PetscCall(FNSetFromOptions(r));
108 2 : np = 2; nq = 3;
109 2 : p[0] = -1.0; p[1] = 0.0;
110 2 : q[0] = 1.0; q[1] = 0.0; q[2] = 1.0;
111 2 : PetscCall(FNRationalSetNumerator(r,np,p));
112 2 : PetscCall(FNRationalSetDenominator(r,nq,q));
113 : /* h(x) */
114 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&h));
115 2 : PetscCall(FNSetType(h,FNCOMBINE));
116 2 : PetscCall(FNSetFromOptions(h));
117 2 : PetscCall(FNCombineSetChildren(h,FN_COMBINE_COMPOSE,r,e));
118 : /* g(x) = 1-x^2 */
119 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&g));
120 2 : PetscCall(FNSetType(g,FNRATIONAL));
121 2 : PetscCall(FNSetFromOptions(g));
122 2 : np = 3;
123 2 : p[0] = -1.0; p[1] = 0.0; p[2] = 1.0;
124 2 : PetscCall(FNRationalSetNumerator(g,np,p));
125 : /* f(x) */
126 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f));
127 2 : PetscCall(FNSetType(f,FNCOMBINE));
128 2 : PetscCall(FNSetFromOptions(f));
129 2 : PetscCall(FNCombineSetChildren(f,FN_COMBINE_MULTIPLY,g,h));
130 :
131 : /* Set up viewer */
132 2 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
133 2 : PetscCall(FNView(f,viewer));
134 2 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
135 :
136 : /* Scalar evaluation */
137 2 : x = 2.2;
138 2 : PetscCall(SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE));
139 2 : PetscCall(FNEvaluateFunction(f,x,&y));
140 2 : PetscCall(FNEvaluateDerivative(f,x,&yp));
141 2 : PetscCall(SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE));
142 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str));
143 2 : PetscCall(SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE));
144 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str));
145 :
146 : /* Test duplication */
147 2 : PetscCall(FNDuplicate(f,PetscObjectComm((PetscObject)f),&fcopy));
148 :
149 : /* Create matrices */
150 2 : if (matcuda) {
151 : #if defined(PETSC_HAVE_CUDA)
152 : PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A));
153 : #endif
154 2 : } else PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A));
155 2 : PetscCall(PetscObjectSetName((PetscObject)A,"A"));
156 :
157 : /* Fill A with a symmetric Toeplitz matrix */
158 2 : PetscCall(MatDenseGetArray(A,&As));
159 22 : for (i=0;i<n;i++) As[i+i*n]=2.0;
160 6 : for (j=1;j<3;j++) {
161 38 : for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
162 : }
163 2 : PetscCall(MatDenseRestoreArray(A,&As));
164 2 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
165 2 : PetscCall(TestMatCombine(fcopy,A,viewer,verbose,inplace));
166 :
167 : /* Repeat with same matrix as non-symmetric */
168 2 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE));
169 2 : PetscCall(TestMatCombine(fcopy,A,viewer,verbose,inplace));
170 :
171 2 : PetscCall(MatDestroy(&A));
172 2 : PetscCall(FNDestroy(&f));
173 2 : PetscCall(FNDestroy(&fcopy));
174 2 : PetscCall(FNDestroy(&g));
175 2 : PetscCall(FNDestroy(&h));
176 2 : PetscCall(FNDestroy(&e));
177 2 : PetscCall(FNDestroy(&r));
178 2 : PetscCall(SlepcFinalize());
179 : return 0;
180 : }
181 :
182 : /*TEST
183 :
184 : testset:
185 : output_file: output/test6_1.out
186 : test:
187 : suffix: 1
188 : test:
189 : suffix: 1_cuda
190 : args: -matcuda
191 : requires: cuda
192 : test:
193 : suffix: 2
194 : args: -inplace
195 : test:
196 : suffix: 2_cuda
197 : args: -inplace -matcuda
198 : requires: cuda
199 :
200 : TEST*/
|