Actual source code: test13.c

slepc-3.17.0 2022-03-31
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 logarithm.\n\n";

 13: #include <slepcfn.h>

 15: /*
 16:    Compute matrix logarithm B = logm(A)
 17:  */
 18: PetscErrorCode TestMatLog(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
 19: {
 20:   PetscBool      set,flg;
 21:   PetscScalar    tau,eta;
 22:   PetscInt       n;
 23:   Mat            F,R;
 24:   Vec            v,f0;
 25:   FN             fnexp;
 26:   PetscReal      nrm;

 29:   MatGetSize(A,&n,NULL);
 30:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
 31:   PetscObjectSetName((PetscObject)F,"F");
 32:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&R);
 33:   PetscObjectSetName((PetscObject)R,"R");
 34:   FNGetScale(fn,&tau,&eta);
 35:   /* compute matrix logarithm */
 36:   if (inplace) {
 37:     MatCopy(A,F,SAME_NONZERO_PATTERN);
 38:     MatIsHermitianKnown(A,&set,&flg);
 39:     if (set && flg) MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE);
 40:     FNEvaluateFunctionMat(fn,F,NULL);
 41:   } else FNEvaluateFunctionMat(fn,A,F);
 42:   if (verbose) {
 43:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 44:     MatView(A,viewer);
 45:     PetscPrintf(PETSC_COMM_WORLD,"Computed logm(A) - - - - - - -\n");
 46:     MatView(F,viewer);
 47:   }
 48:   /* check error ||expm(F)-A||_F */
 49:   FNCreate(PETSC_COMM_WORLD,&fnexp);
 50:   FNSetType(fnexp,FNEXP);
 51:   MatCopy(F,R,SAME_NONZERO_PATTERN);
 52:   if (eta!=1.0) MatScale(R,1.0/eta);
 53:   FNEvaluateFunctionMat(fnexp,R,NULL);
 54:   FNDestroy(&fnexp);
 55:   MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
 56:   MatNorm(R,NORM_FROBENIUS,&nrm);
 57:   if (nrm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"||expm(F)-A||_F < 100*eps\n");
 58:   else PetscPrintf(PETSC_COMM_WORLD,"||expm(F)-A||_F = %g\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:   MatDestroy(&R);
 68:   VecDestroy(&v);
 69:   VecDestroy(&f0);
 70:   PetscFunctionReturn(0);
 71: }

 73: int main(int argc,char **argv)
 74: {
 75:   FN             fn;
 76:   Mat            A;
 77:   PetscInt       i,j,n=10;
 78:   PetscScalar    *As;
 79:   PetscViewer    viewer;
 80:   PetscBool      verbose,inplace,random,triang;

 82:   SlepcInitialize(&argc,&argv,(char*)0,help);
 83:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 84:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 85:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
 86:   PetscOptionsHasName(NULL,NULL,"-random",&random);
 87:   PetscOptionsHasName(NULL,NULL,"-triang",&triang);
 88:   PetscPrintf(PETSC_COMM_WORLD,"Matrix logarithm, n=%" PetscInt_FMT ".\n",n);

 90:   /* Create logarithm function object */
 91:   FNCreate(PETSC_COMM_WORLD,&fn);
 92:   FNSetType(fn,FNLOG);
 93:   FNSetFromOptions(fn);

 95:   /* Set up viewer */
 96:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 97:   FNView(fn,viewer);
 98:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

100:   /* Create matrices */
101:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
102:   PetscObjectSetName((PetscObject)A,"A");

104:   if (random) MatSetRandom(A,NULL);
105:   else {
106:     /* Fill A with a non-symmetric Toeplitz matrix */
107:     MatDenseGetArray(A,&As);
108:     for (i=0;i<n;i++) As[i+i*n]=2.0;
109:     for (j=1;j<3;j++) {
110:       for (i=0;i<n-j;i++) {
111:         As[i+(i+j)*n]=1.0;
112:         if (!triang) As[(i+j)+i*n]=-1.0;
113:       }
114:     }
115:     As[(n-1)*n] = -5.0;
116:     As[0] = 2.01;
117:     MatDenseRestoreArray(A,&As);
118:   }
119:   TestMatLog(fn,A,viewer,verbose,inplace);

121:   MatDestroy(&A);
122:   FNDestroy(&fn);
123:   SlepcFinalize();
124:   return 0;
125: }

127: /*TEST

129:    testset:
130:       filter: grep -v "computing matrix functions"
131:       output_file: output/test13_1.out
132:       test:
133:          suffix: 1
134:          args: -fn_scale .04,2 -n 75
135:          requires: c99_complex !__float128
136:       test:
137:          suffix: 1_triang
138:          args: -fn_scale .04,2 -n 75 -triang
139:          requires: c99_complex !__float128
140:       test:
141:          suffix: 1_random
142:          args: -fn_scale .02,2 -n 75 -random
143:          requires: complex
144:          filter_output: sed -e 's/04/02/'

146: TEST*/