Actual source code: test5.c

slepc-3.21.2 2024-09-25
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 rational function.\n\n";

 13: #include <slepcfn.h>

 15: /*
 16:    Compute matrix rational function B = q(A)\p(A)
 17:  */
 18: PetscErrorCode TestMatRational(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
 19: {
 20:   PetscBool      set,flg;
 21:   PetscInt       n;
 22:   Mat            F,Acopy;
 23:   Vec            v,f0;
 24:   PetscReal      nrm;

 26:   PetscFunctionBeginUser;
 27:   PetscCall(MatGetSize(A,&n,NULL));
 28:   PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&F));
 29:   PetscCall(PetscObjectSetName((PetscObject)F,"F"));
 30:   /* compute matrix function */
 31:   if (inplace) {
 32:     PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN));
 33:     PetscCall(MatIsHermitianKnown(A,&set,&flg));
 34:     if (set && flg) PetscCall(MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE));
 35:     PetscCall(FNEvaluateFunctionMat(fn,F,NULL));
 36:   } else {
 37:     PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Acopy));
 38:     PetscCall(FNEvaluateFunctionMat(fn,A,F));
 39:     /* check that A has not been modified */
 40:     PetscCall(MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN));
 41:     PetscCall(MatNorm(Acopy,NORM_1,&nrm));
 42:     if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm));
 43:     PetscCall(MatDestroy(&Acopy));
 44:   }
 45:   if (verbose) {
 46:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n"));
 47:     PetscCall(MatView(A,viewer));
 48:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n"));
 49:     PetscCall(MatView(F,viewer));
 50:   }
 51:   /* print matrix norm for checking */
 52:   PetscCall(MatNorm(F,NORM_1,&nrm));
 53:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm));
 54:   /* check FNEvaluateFunctionMatVec() */
 55:   PetscCall(MatCreateVecs(A,&v,&f0));
 56:   PetscCall(MatGetColumnVector(F,f0,0));
 57:   PetscCall(FNEvaluateFunctionMatVec(fn,A,v));
 58:   PetscCall(VecAXPY(v,-1.0,f0));
 59:   PetscCall(VecNorm(v,NORM_2,&nrm));
 60:   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));
 61:   PetscCall(MatDestroy(&F));
 62:   PetscCall(VecDestroy(&v));
 63:   PetscCall(VecDestroy(&f0));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: int main(int argc,char **argv)
 68: {
 69:   FN             fn;
 70:   Mat            A=NULL;
 71:   PetscInt       i,j,n=10,np,nq;
 72:   PetscScalar    *As,p[10],q[10];
 73:   PetscViewer    viewer;
 74:   PetscBool      verbose,inplace,matcuda;

 76:   PetscFunctionBeginUser;
 77:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 78:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 79:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 80:   PetscCall(PetscOptionsHasName(NULL,NULL,"-inplace",&inplace));
 81:   PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda));
 82:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix rational function, n=%" PetscInt_FMT ".\n",n));

 84:   /* Create rational function r(x)=p(x)/q(x) */
 85:   PetscCall(FNCreate(PETSC_COMM_WORLD,&fn));
 86:   PetscCall(FNSetType(fn,FNRATIONAL));
 87:   np = 2; nq = 3;
 88:   p[0] = -3.1; p[1] = 1.1;
 89:   q[0] = 1.0; q[1] = -2.0; q[2] = 3.5;
 90:   PetscCall(FNRationalSetNumerator(fn,np,p));
 91:   PetscCall(FNRationalSetDenominator(fn,nq,q));
 92:   PetscCall(FNSetFromOptions(fn));

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

 99:   /* Create matrices */
100:   if (matcuda) {
101: #if defined(PETSC_HAVE_CUDA)
102:     PetscCall(MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A));
103: #endif
104:   } else PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A));
105:   PetscCall(PetscObjectSetName((PetscObject)A,"A"));

107:   /* Fill A with a symmetric Toeplitz matrix */
108:   PetscCall(MatDenseGetArray(A,&As));
109:   for (i=0;i<n;i++) As[i+i*n]=2.0;
110:   for (j=1;j<3;j++) {
111:     for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
112:   }
113:   PetscCall(MatDenseRestoreArray(A,&As));
114:   PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
115:   PetscCall(TestMatRational(fn,A,viewer,verbose,inplace));

117:   /* Repeat with same matrix as non-symmetric */
118:   PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE));
119:   PetscCall(TestMatRational(fn,A,viewer,verbose,inplace));

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

127: /*TEST

129:    testset:
130:       output_file: output/test5_1.out
131:       requires: !single
132:       test:
133:          suffix: 1
134:       test:
135:          suffix: 1_cuda
136:          args: -matcuda
137:          requires: cuda
138:       test:
139:          suffix: 2
140:          args: -inplace
141:       test:
142:          suffix: 2_cuda
143:          args: -inplace -matcuda
144:          requires: cuda

146: TEST*/