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 : static char help[] = "Test matrix square root.\n\n";
12 :
13 : #include <slepcfn.h>
14 :
15 : /*
16 : Compute matrix square root B = sqrtm(A)
17 : Check result as norm(B*B-A)
18 : */
19 42 : PetscErrorCode TestMatSqrt(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
20 : {
21 42 : PetscScalar tau,eta;
22 42 : PetscReal nrm;
23 42 : PetscBool set,flg;
24 42 : PetscInt n;
25 42 : Mat S,R,Acopy;
26 42 : Vec v,f0;
27 :
28 42 : PetscFunctionBeginUser;
29 42 : PetscCall(MatGetSize(A,&n,NULL));
30 42 : PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&S));
31 42 : PetscCall(PetscObjectSetName((PetscObject)S,"S"));
32 42 : PetscCall(FNGetScale(fn,&tau,&eta));
33 : /* compute square root */
34 42 : if (inplace) {
35 21 : PetscCall(MatCopy(A,S,SAME_NONZERO_PATTERN));
36 21 : PetscCall(MatIsHermitianKnown(A,&set,&flg));
37 21 : if (set && flg) PetscCall(MatSetOption(S,MAT_HERMITIAN,PETSC_TRUE));
38 21 : PetscCall(FNEvaluateFunctionMat(fn,S,NULL));
39 : } else {
40 21 : PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Acopy));
41 21 : PetscCall(FNEvaluateFunctionMat(fn,A,S));
42 : /* check that A has not been modified */
43 21 : PetscCall(MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN));
44 21 : PetscCall(MatNorm(Acopy,NORM_1,&nrm));
45 21 : if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm));
46 21 : PetscCall(MatDestroy(&Acopy));
47 : }
48 42 : if (verbose) {
49 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n"));
50 0 : PetscCall(MatView(A,viewer));
51 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed sqrtm(A) - - - - - - -\n"));
52 0 : PetscCall(MatView(S,viewer));
53 : }
54 : /* check error ||S*S-A||_F */
55 42 : PetscCall(MatMatMult(S,S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R));
56 42 : if (eta!=1.0) PetscCall(MatScale(R,1.0/(eta*eta)));
57 42 : PetscCall(MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN));
58 42 : PetscCall(MatNorm(R,NORM_FROBENIUS,&nrm));
59 42 : if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n"));
60 4 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm));
61 : /* check FNEvaluateFunctionMatVec() */
62 42 : PetscCall(MatCreateVecs(A,&v,&f0));
63 42 : PetscCall(MatGetColumnVector(S,f0,0));
64 42 : PetscCall(FNEvaluateFunctionMatVec(fn,A,v));
65 42 : PetscCall(VecAXPY(v,-1.0,f0));
66 42 : PetscCall(VecNorm(v,NORM_2,&nrm));
67 42 : 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));
68 42 : PetscCall(MatDestroy(&S));
69 42 : PetscCall(MatDestroy(&R));
70 42 : PetscCall(VecDestroy(&v));
71 42 : PetscCall(VecDestroy(&f0));
72 42 : PetscFunctionReturn(PETSC_SUCCESS);
73 : }
74 :
75 14 : int main(int argc,char **argv)
76 : {
77 14 : FN fn;
78 14 : Mat A=NULL;
79 14 : PetscInt i,j,n=10;
80 14 : PetscScalar *As;
81 14 : PetscViewer viewer;
82 14 : PetscBool verbose,inplace,matcuda;
83 14 : PetscRandom myrand;
84 14 : PetscReal v;
85 :
86 14 : PetscFunctionBeginUser;
87 14 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
88 14 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
89 14 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
90 14 : PetscCall(PetscOptionsHasName(NULL,NULL,"-inplace",&inplace));
91 14 : PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda));
92 14 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%" PetscInt_FMT ".\n",n));
93 :
94 : /* Create function object */
95 14 : PetscCall(FNCreate(PETSC_COMM_WORLD,&fn));
96 14 : PetscCall(FNSetType(fn,FNSQRT));
97 14 : PetscCall(FNSetFromOptions(fn));
98 :
99 : /* Set up viewer */
100 14 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
101 14 : PetscCall(FNView(fn,viewer));
102 14 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
103 :
104 : /* Create matrix */
105 14 : if (matcuda) {
106 : #if defined(PETSC_HAVE_CUDA)
107 : PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A));
108 : #endif
109 14 : } else PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A));
110 14 : PetscCall(PetscObjectSetName((PetscObject)A,"A"));
111 :
112 : /* Compute square root of a symmetric matrix A */
113 14 : PetscCall(MatDenseGetArray(A,&As));
114 1414 : for (i=0;i<n;i++) As[i+i*n]=2.5;
115 42 : for (j=1;j<3;j++) {
116 2786 : for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
117 : }
118 14 : PetscCall(MatDenseRestoreArray(A,&As));
119 14 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
120 14 : PetscCall(TestMatSqrt(fn,A,viewer,verbose,inplace));
121 :
122 : /* Repeat with upper triangular A */
123 14 : PetscCall(MatDenseGetArray(A,&As));
124 42 : for (j=1;j<3;j++) {
125 2786 : for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
126 : }
127 14 : PetscCall(MatDenseRestoreArray(A,&As));
128 14 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE));
129 14 : PetscCall(TestMatSqrt(fn,A,viewer,verbose,inplace));
130 :
131 : /* Repeat with non-symmetic A */
132 14 : PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&myrand));
133 14 : PetscCall(PetscRandomSetFromOptions(myrand));
134 14 : PetscCall(PetscRandomSetInterval(myrand,0.0,1.0));
135 14 : PetscCall(MatDenseGetArray(A,&As));
136 42 : for (j=1;j<3;j++) {
137 2786 : for (i=0;i<n-j;i++) {
138 2758 : PetscCall(PetscRandomGetValueReal(myrand,&v));
139 2758 : As[(i+j)+i*n]=v;
140 : }
141 : }
142 14 : PetscCall(MatDenseRestoreArray(A,&As));
143 14 : PetscCall(PetscRandomDestroy(&myrand));
144 14 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE));
145 14 : PetscCall(TestMatSqrt(fn,A,viewer,verbose,inplace));
146 :
147 14 : PetscCall(MatDestroy(&A));
148 14 : PetscCall(FNDestroy(&fn));
149 14 : PetscCall(SlepcFinalize());
150 : return 0;
151 : }
152 :
153 : /*TEST
154 :
155 : testset:
156 : args: -fn_scale .05,2 -n 100
157 : filter: grep -v "computing matrix functions"
158 : output_file: output/test7_1.out
159 : requires: !__float128
160 : timeoutfactor: 2
161 : test:
162 : suffix: 1
163 : args: -fn_method {{0 1 2}}
164 : test:
165 : suffix: 1_sadeghi
166 : args: -fn_method 3
167 : requires: !single
168 : test:
169 : suffix: 1_cuda
170 : args: -fn_method 2 -matcuda
171 : requires: cuda !single
172 : test:
173 : suffix: 1_magma
174 : args: -fn_method {{1 3}} -matcuda
175 : requires: cuda magma !single
176 : test:
177 : suffix: 2
178 : args: -inplace -fn_method {{0 1 2}}
179 : test:
180 : suffix: 2_sadeghi
181 : args: -inplace -fn_method 3
182 : requires: !single
183 : test:
184 : suffix: 2_cuda
185 : args: -inplace -fn_method 2 -matcuda
186 : requires: cuda !single
187 : test:
188 : suffix: 2_magma
189 : args: -inplace -fn_method {{1 3}} -matcuda
190 : requires: cuda magma !single
191 :
192 : testset:
193 : nsize: 3
194 : args: -fn_scale .05,2 -n 100 -fn_parallel synchronized
195 : filter: grep -v "computing matrix functions" | grep -v "SYNCHRONIZED" | sed -e "s/3 MPI processes/1 MPI process/g"
196 : requires: !__float128
197 : output_file: output/test7_1.out
198 : test:
199 : suffix: 3
200 : test:
201 : suffix: 3_inplace
202 : args: -inplace
203 :
204 : TEST*/
|