Actual source code: test7.c

slepc-3.18.3 2023-03-24
Report Typos and Errors
  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;

 29:   MatGetSize(A,&n,NULL);
 30:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&S);
 31:   PetscObjectSetName((PetscObject)S,"S");
 32:   FNGetScale(fn,&tau,&eta);
 33:   /* compute square root */
 34:   if (inplace) {
 35:     MatCopy(A,S,SAME_NONZERO_PATTERN);
 36:     MatIsHermitianKnown(A,&set,&flg);
 37:     if (set && flg) MatSetOption(S,MAT_HERMITIAN,PETSC_TRUE);
 38:     FNEvaluateFunctionMat(fn,S,NULL);
 39:   } else {
 40:     MatDuplicate(A,MAT_COPY_VALUES,&Acopy);
 41:     FNEvaluateFunctionMat(fn,A,S);
 42:     /* check that A has not been modified */
 43:     MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN);
 44:     MatNorm(Acopy,NORM_1,&nrm);
 45:     if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm);
 46:     MatDestroy(&Acopy);
 47:   }
 48:   if (verbose) {
 49:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 50:     MatView(A,viewer);
 51:     PetscPrintf(PETSC_COMM_WORLD,"Computed sqrtm(A) - - - - - - -\n");
 52:     MatView(S,viewer);
 53:   }
 54:   /* check error ||S*S-A||_F */
 55:   MatMatMult(S,S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
 56:   if (eta!=1.0) MatScale(R,1.0/(eta*eta));
 57:   MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
 58:   MatNorm(R,NORM_FROBENIUS,&nrm);
 59:   if (nrm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n");
 60:   else PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm);
 61:   /* check FNEvaluateFunctionMatVec() */
 62:   MatCreateVecs(A,&v,&f0);
 63:   MatGetColumnVector(S,f0,0);
 64:   FNEvaluateFunctionMatVec(fn,A,v);
 65:   VecAXPY(v,-1.0,f0);
 66:   VecNorm(v,NORM_2,&nrm);
 67:   if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 68:   MatDestroy(&S);
 69:   MatDestroy(&R);
 70:   VecDestroy(&v);
 71:   VecDestroy(&f0);
 72:   return 0;
 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;

 87:   SlepcInitialize(&argc,&argv,(char*)0,help);
 88:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 89:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 90:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
 91:   PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda);
 92:   PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%" PetscInt_FMT ".\n",n);

 94:   /* Create function object */
 95:   FNCreate(PETSC_COMM_WORLD,&fn);
 96:   FNSetType(fn,FNSQRT);
 97:   FNSetFromOptions(fn);

 99:   /* Set up viewer */
100:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
101:   FNView(fn,viewer);
102:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

104:   /* Create matrix */
105:   if (matcuda) {
106: #if defined(PETSC_HAVE_CUDA)
107:     MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A);
108: #endif
109:   } else MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
110:   PetscObjectSetName((PetscObject)A,"A");

112:   /* Compute square root of a symmetric matrix A */
113:   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:   MatDenseRestoreArray(A,&As);
119:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
120:   TestMatSqrt(fn,A,viewer,verbose,inplace);

122:   /* Repeat with upper triangular A */
123:   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:   MatDenseRestoreArray(A,&As);
128:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
129:   TestMatSqrt(fn,A,viewer,verbose,inplace);

131:   /* Repeat with non-symmetic A */
132:   PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
133:   PetscRandomSetFromOptions(myrand);
134:   PetscRandomSetInterval(myrand,0.0,1.0);
135:   MatDenseGetArray(A,&As);
136:   for (j=1;j<3;j++) {
137:     for (i=0;i<n-j;i++) {
138:       PetscRandomGetValueReal(myrand,&v);
139:       As[(i+j)+i*n]=v;
140:     }
141:   }
142:   MatDenseRestoreArray(A,&As);
143:   PetscRandomDestroy(&myrand);
144:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
145:   TestMatSqrt(fn,A,viewer,verbose,inplace);

147:   MatDestroy(&A);
148:   FNDestroy(&fn);
149:   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*/