Actual source code: test6.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: */
 10: /*
 11:    Define the function

 13:         f(x) = (1-x^2) exp(-x/(1+x^2))

 15:    with the following tree:

 17:             f(x)                  f(x)              (combined by product)
 18:            /    \                 g(x) = 1-x^2      (polynomial)
 19:         g(x)    h(x)              h(x)              (combined by composition)
 20:                /    \             r(x) = -x/(1+x^2) (rational)
 21:              r(x)   e(x)          e(x) = exp(x)     (exponential)
 22: */

 24: static char help[] = "Test combined function.\n\n";

 26: #include <slepcfn.h>

 28: /*
 29:    Compute matrix function B = (I-A^2) exp(-(I+A^2)\A)
 30:  */
 31: PetscErrorCode TestMatCombine(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
 32: {
 33:   PetscBool      set,flg;
 34:   PetscInt       n;
 35:   Mat            F;
 36:   Vec            v,f0;
 37:   PetscReal      nrm;

 40:   MatGetSize(A,&n,NULL);
 41:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
 42:   PetscObjectSetName((PetscObject)F,"F");
 43:   /* compute matrix function */
 44:   if (inplace) {
 45:     MatCopy(A,F,SAME_NONZERO_PATTERN);
 46:     MatIsHermitianKnown(A,&set,&flg);
 47:     if (set && flg) MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE);
 48:     FNEvaluateFunctionMat(fn,F,NULL);
 49:   } else FNEvaluateFunctionMat(fn,A,F);
 50:   if (verbose) {
 51:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 52:     MatView(A,viewer);
 53:     PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
 54:     MatView(F,viewer);
 55:   }
 56:   /* print matrix norm for checking */
 57:   MatNorm(F,NORM_1,&nrm);
 58:   PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %6.3f\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:   VecDestroy(&v);
 68:   VecDestroy(&f0);
 69:   PetscFunctionReturn(0);
 70: }

 72: int main(int argc,char **argv)
 73: {
 74:   FN             f,g,h,e,r,fcopy;
 75:   Mat            A;
 76:   PetscInt       i,j,n=10,np,nq;
 77:   PetscScalar    x,y,yp,*As,p[10],q[10];
 78:   char           strx[50],str[50];
 79:   PetscViewer    viewer;
 80:   PetscBool      verbose,inplace;

 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:   PetscPrintf(PETSC_COMM_WORLD,"Combined function, n=%" PetscInt_FMT ".\n",n);

 88:   /* Create function */

 90:   /* e(x) = exp(x) */
 91:   FNCreate(PETSC_COMM_WORLD,&e);
 92:   FNSetType(e,FNEXP);
 93:   FNSetFromOptions(e);
 94:   /* r(x) = x/(1+x^2) */
 95:   FNCreate(PETSC_COMM_WORLD,&r);
 96:   FNSetType(r,FNRATIONAL);
 97:   FNSetFromOptions(r);
 98:   np = 2; nq = 3;
 99:   p[0] = -1.0; p[1] = 0.0;
100:   q[0] = 1.0; q[1] = 0.0; q[2] = 1.0;
101:   FNRationalSetNumerator(r,np,p);
102:   FNRationalSetDenominator(r,nq,q);
103:   /* h(x) */
104:   FNCreate(PETSC_COMM_WORLD,&h);
105:   FNSetType(h,FNCOMBINE);
106:   FNSetFromOptions(h);
107:   FNCombineSetChildren(h,FN_COMBINE_COMPOSE,r,e);
108:   /* g(x) = 1-x^2 */
109:   FNCreate(PETSC_COMM_WORLD,&g);
110:   FNSetType(g,FNRATIONAL);
111:   FNSetFromOptions(g);
112:   np = 3;
113:   p[0] = -1.0; p[1] = 0.0; p[2] = 1.0;
114:   FNRationalSetNumerator(g,np,p);
115:   /* f(x) */
116:   FNCreate(PETSC_COMM_WORLD,&f);
117:   FNSetType(f,FNCOMBINE);
118:   FNSetFromOptions(f);
119:   FNCombineSetChildren(f,FN_COMBINE_MULTIPLY,g,h);

121:   /* Set up viewer */
122:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
123:   FNView(f,viewer);
124:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

126:   /* Scalar evaluation */
127:   x = 2.2;
128:   SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
129:   FNEvaluateFunction(f,x,&y);
130:   FNEvaluateDerivative(f,x,&yp);
131:   SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
132:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
133:   SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
134:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

136:   /* Test duplication */
137:   FNDuplicate(f,PetscObjectComm((PetscObject)f),&fcopy);

139:   /* Create matrices */
140:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
141:   PetscObjectSetName((PetscObject)A,"A");

143:   /* Fill A with a symmetric Toeplitz matrix */
144:   MatDenseGetArray(A,&As);
145:   for (i=0;i<n;i++) As[i+i*n]=2.0;
146:   for (j=1;j<3;j++) {
147:     for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
148:   }
149:   MatDenseRestoreArray(A,&As);
150:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
151:   TestMatCombine(fcopy,A,viewer,verbose,inplace);

153:   /* Repeat with same matrix as non-symmetric */
154:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
155:   TestMatCombine(fcopy,A,viewer,verbose,inplace);

157:   MatDestroy(&A);
158:   FNDestroy(&f);
159:   FNDestroy(&fcopy);
160:   FNDestroy(&g);
161:   FNDestroy(&h);
162:   FNDestroy(&e);
163:   FNDestroy(&r);
164:   SlepcFinalize();
165:   return 0;
166: }

168: /*TEST

170:    test:
171:       suffix: 1
172:       nsize: 1

174:    test:
175:       suffix: 2
176:       nsize: 1
177:       args: -inplace
178:       output_file: output/test6_1.out

180: TEST*/