Actual source code: test7.c
slepc-3.21.1 2024-04-26
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: */
11: static char help[] = "Test matrix square root.\n\n";
13: #include <slepcfn.h>
15: /*
16: Compute matrix square root B = sqrtm(A)
17: Check result as norm(B*B-A)
18: */
19: PetscErrorCode TestMatSqrt(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
20: {
21: PetscScalar tau,eta;
22: PetscReal nrm;
23: PetscBool set,flg;
24: PetscInt n;
25: Mat S,R,Acopy;
26: Vec v,f0;
28: PetscFunctionBeginUser;
29: PetscCall(MatGetSize(A,&n,NULL));
30: PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&S));
31: PetscCall(PetscObjectSetName((PetscObject)S,"S"));
32: PetscCall(FNGetScale(fn,&tau,&eta));
33: /* compute square root */
34: if (inplace) {
35: PetscCall(MatCopy(A,S,SAME_NONZERO_PATTERN));
36: PetscCall(MatIsHermitianKnown(A,&set,&flg));
37: if (set && flg) PetscCall(MatSetOption(S,MAT_HERMITIAN,PETSC_TRUE));
38: PetscCall(FNEvaluateFunctionMat(fn,S,NULL));
39: } else {
40: PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Acopy));
41: PetscCall(FNEvaluateFunctionMat(fn,A,S));
42: /* check that A has not been modified */
43: PetscCall(MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN));
44: PetscCall(MatNorm(Acopy,NORM_1,&nrm));
45: if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm));
46: PetscCall(MatDestroy(&Acopy));
47: }
48: if (verbose) {
49: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n"));
50: PetscCall(MatView(A,viewer));
51: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed sqrtm(A) - - - - - - -\n"));
52: PetscCall(MatView(S,viewer));
53: }
54: /* check error ||S*S-A||_F */
55: PetscCall(MatMatMult(S,S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R));
56: if (eta!=1.0) PetscCall(MatScale(R,1.0/(eta*eta)));
57: PetscCall(MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN));
58: PetscCall(MatNorm(R,NORM_FROBENIUS,&nrm));
59: if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n"));
60: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm));
61: /* check FNEvaluateFunctionMatVec() */
62: PetscCall(MatCreateVecs(A,&v,&f0));
63: PetscCall(MatGetColumnVector(S,f0,0));
64: PetscCall(FNEvaluateFunctionMatVec(fn,A,v));
65: PetscCall(VecAXPY(v,-1.0,f0));
66: PetscCall(VecNorm(v,NORM_2,&nrm));
67: 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: PetscCall(MatDestroy(&S));
69: PetscCall(MatDestroy(&R));
70: PetscCall(VecDestroy(&v));
71: PetscCall(VecDestroy(&f0));
72: PetscFunctionReturn(PETSC_SUCCESS);
73: }
75: int main(int argc,char **argv)
76: {
77: FN fn;
78: Mat A=NULL;
79: PetscInt i,j,n=10;
80: PetscScalar *As;
81: PetscViewer viewer;
82: PetscBool verbose,inplace,matcuda;
83: PetscRandom myrand;
84: PetscReal v;
86: PetscFunctionBeginUser;
87: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
88: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
89: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
90: PetscCall(PetscOptionsHasName(NULL,NULL,"-inplace",&inplace));
91: PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda));
92: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%" PetscInt_FMT ".\n",n));
94: /* Create function object */
95: PetscCall(FNCreate(PETSC_COMM_WORLD,&fn));
96: PetscCall(FNSetType(fn,FNSQRT));
97: PetscCall(FNSetFromOptions(fn));
99: /* Set up viewer */
100: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
101: PetscCall(FNView(fn,viewer));
102: if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
104: /* Create matrix */
105: if (matcuda) {
106: #if defined(PETSC_HAVE_CUDA)
107: PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A));
108: #endif
109: } else PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A));
110: PetscCall(PetscObjectSetName((PetscObject)A,"A"));
112: /* Compute square root of a symmetric matrix A */
113: PetscCall(MatDenseGetArray(A,&As));
114: for (i=0;i<n;i++) As[i+i*n]=2.5;
115: for (j=1;j<3;j++) {
116: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
117: }
118: PetscCall(MatDenseRestoreArray(A,&As));
119: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
120: PetscCall(TestMatSqrt(fn,A,viewer,verbose,inplace));
122: /* Repeat with upper triangular A */
123: PetscCall(MatDenseGetArray(A,&As));
124: for (j=1;j<3;j++) {
125: for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
126: }
127: PetscCall(MatDenseRestoreArray(A,&As));
128: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE));
129: PetscCall(TestMatSqrt(fn,A,viewer,verbose,inplace));
131: /* Repeat with non-symmetic A */
132: PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&myrand));
133: PetscCall(PetscRandomSetFromOptions(myrand));
134: PetscCall(PetscRandomSetInterval(myrand,0.0,1.0));
135: PetscCall(MatDenseGetArray(A,&As));
136: for (j=1;j<3;j++) {
137: for (i=0;i<n-j;i++) {
138: PetscCall(PetscRandomGetValueReal(myrand,&v));
139: As[(i+j)+i*n]=v;
140: }
141: }
142: PetscCall(MatDenseRestoreArray(A,&As));
143: PetscCall(PetscRandomDestroy(&myrand));
144: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE));
145: PetscCall(TestMatSqrt(fn,A,viewer,verbose,inplace));
147: PetscCall(MatDestroy(&A));
148: PetscCall(FNDestroy(&fn));
149: PetscCall(SlepcFinalize());
150: return 0;
151: }
153: /*TEST
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
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
204: TEST*/