Actual source code: test6.c

slepc-3.18.1 2022-11-02
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,Acopy;
 36:   Vec            v,f0;
 37:   PetscReal      nrm;

 40:   MatGetSize(A,&n,NULL);
 41:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&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 {
 50:     MatDuplicate(A,MAT_COPY_VALUES,&Acopy);
 51:     FNEvaluateFunctionMat(fn,A,F);
 52:     /* check that A has not been modified */
 53:     MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN);
 54:     MatNorm(Acopy,NORM_1,&nrm);
 55:     if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm);
 56:     MatDestroy(&Acopy);
 57:   }
 58:   if (verbose) {
 59:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 60:     MatView(A,viewer);
 61:     PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
 62:     MatView(F,viewer);
 63:   }
 64:   /* print matrix norm for checking */
 65:   MatNorm(F,NORM_1,&nrm);
 66:   PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %6.3f\n",(double)nrm);
 67:   /* check FNEvaluateFunctionMatVec() */
 68:   MatCreateVecs(A,&v,&f0);
 69:   MatGetColumnVector(F,f0,0);
 70:   FNEvaluateFunctionMatVec(fn,A,v);
 71:   VecAXPY(v,-1.0,f0);
 72:   VecNorm(v,NORM_2,&nrm);
 73:   if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 74:   MatDestroy(&F);
 75:   VecDestroy(&v);
 76:   VecDestroy(&f0);
 77:   return 0;
 78: }

 80: int main(int argc,char **argv)
 81: {
 82:   FN             f,g,h,e,r,fcopy;
 83:   Mat            A=NULL;
 84:   PetscInt       i,j,n=10,np,nq;
 85:   PetscScalar    x,y,yp,*As,p[10],q[10];
 86:   char           strx[50],str[50];
 87:   PetscViewer    viewer;
 88:   PetscBool      verbose,inplace,matcuda;

 91:   SlepcInitialize(&argc,&argv,(char*)0,help);
 92:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 93:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 94:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
 95:   PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda);
 96:   PetscPrintf(PETSC_COMM_WORLD,"Combined function, n=%" PetscInt_FMT ".\n",n);

 98:   /* Create function */

100:   /* e(x) = exp(x) */
101:   FNCreate(PETSC_COMM_WORLD,&e);
102:   FNSetType(e,FNEXP);
103:   FNSetFromOptions(e);
104:   /* r(x) = x/(1+x^2) */
105:   FNCreate(PETSC_COMM_WORLD,&r);
106:   FNSetType(r,FNRATIONAL);
107:   FNSetFromOptions(r);
108:   np = 2; nq = 3;
109:   p[0] = -1.0; p[1] = 0.0;
110:   q[0] = 1.0; q[1] = 0.0; q[2] = 1.0;
111:   FNRationalSetNumerator(r,np,p);
112:   FNRationalSetDenominator(r,nq,q);
113:   /* h(x) */
114:   FNCreate(PETSC_COMM_WORLD,&h);
115:   FNSetType(h,FNCOMBINE);
116:   FNSetFromOptions(h);
117:   FNCombineSetChildren(h,FN_COMBINE_COMPOSE,r,e);
118:   /* g(x) = 1-x^2 */
119:   FNCreate(PETSC_COMM_WORLD,&g);
120:   FNSetType(g,FNRATIONAL);
121:   FNSetFromOptions(g);
122:   np = 3;
123:   p[0] = -1.0; p[1] = 0.0; p[2] = 1.0;
124:   FNRationalSetNumerator(g,np,p);
125:   /* f(x) */
126:   FNCreate(PETSC_COMM_WORLD,&f);
127:   FNSetType(f,FNCOMBINE);
128:   FNSetFromOptions(f);
129:   FNCombineSetChildren(f,FN_COMBINE_MULTIPLY,g,h);

131:   /* Set up viewer */
132:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
133:   FNView(f,viewer);
134:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

136:   /* Scalar evaluation */
137:   x = 2.2;
138:   SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
139:   FNEvaluateFunction(f,x,&y);
140:   FNEvaluateDerivative(f,x,&yp);
141:   SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
142:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
143:   SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
144:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

146:   /* Test duplication */
147:   FNDuplicate(f,PetscObjectComm((PetscObject)f),&fcopy);

149:   /* Create matrices */
150:   if (matcuda) {
151: #if defined(PETSC_HAVE_CUDA)
152:     MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A);
153: #endif
154:   } else MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
155:   PetscObjectSetName((PetscObject)A,"A");

157:   /* Fill A with a symmetric Toeplitz matrix */
158:   MatDenseGetArray(A,&As);
159:   for (i=0;i<n;i++) As[i+i*n]=2.0;
160:   for (j=1;j<3;j++) {
161:     for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
162:   }
163:   MatDenseRestoreArray(A,&As);
164:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
165:   TestMatCombine(fcopy,A,viewer,verbose,inplace);

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

171:   MatDestroy(&A);
172:   FNDestroy(&f);
173:   FNDestroy(&fcopy);
174:   FNDestroy(&g);
175:   FNDestroy(&h);
176:   FNDestroy(&e);
177:   FNDestroy(&r);
178:   SlepcFinalize();
179:   return 0;
180: }

182: /*TEST

184:    testset:
185:       output_file: output/test6_1.out
186:       test:
187:          suffix: 1
188:       test:
189:          suffix: 1_cuda
190:          args: -matcuda
191:          requires: cuda
192:       test:
193:          suffix: 2
194:          args: -inplace
195:       test:
196:          suffix: 2_cuda
197:          args: -inplace -matcuda
198:          requires: cuda

200: TEST*/