Actual source code: test1.c

slepc-3.17.2 2022-08-09
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 LME interface functions, based on ex32.c.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepclme.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat                  A,B,C,C1,D;
 21:   LME                  lme;
 22:   PetscReal            tol,errest,error;
 23:   PetscScalar          *u;
 24:   PetscInt             N,n=10,m,Istart,Iend,II,maxit,ncv,i,j;
 25:   PetscBool            flg,testprefix=PETSC_FALSE,viewmatrices=PETSC_FALSE;
 26:   const char           *prefix;
 27:   LMEType              type;
 28:   LMEProblemType       ptype;
 29:   PetscViewerAndFormat *vf;

 31:   SlepcInitialize(&argc,&argv,(char*)0,help);

 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 34:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flg);
 35:   if (!flg) m=n;
 36:   N = n*m;
 37:   PetscPrintf(PETSC_COMM_WORLD,"\nLyapunov equation, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
 38:   PetscOptionsGetBool(NULL,NULL,"-test_prefix",&testprefix,NULL);
 39:   PetscOptionsGetBool(NULL,NULL,"-view_matrices",&viewmatrices,NULL);

 41:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 42:                        Create the 2-D Laplacian, A
 43:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 45:   MatCreate(PETSC_COMM_WORLD,&A);
 46:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 47:   MatSetFromOptions(A);
 48:   MatSetUp(A);
 49:   MatGetOwnershipRange(A,&Istart,&Iend);
 50:   for (II=Istart;II<Iend;II++) {
 51:     i = II/n; j = II-i*n;
 52:     if (i>0) MatSetValue(A,II,II-n,1.0,INSERT_VALUES);
 53:     if (i<m-1) MatSetValue(A,II,II+n,1.0,INSERT_VALUES);
 54:     if (j>0) MatSetValue(A,II,II-1,1.0,INSERT_VALUES);
 55:     if (j<n-1) MatSetValue(A,II,II+1,1.0,INSERT_VALUES);
 56:     MatSetValue(A,II,II,-4.0,INSERT_VALUES);
 57:   }
 58:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 62:        Create a low-rank Mat to store the right-hand side C = C1*C1'
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   MatCreate(PETSC_COMM_WORLD,&C1);
 66:   MatSetSizes(C1,PETSC_DECIDE,PETSC_DECIDE,N,2);
 67:   MatSetType(C1,MATDENSE);
 68:   MatSetUp(C1);
 69:   MatGetOwnershipRange(C1,&Istart,&Iend);
 70:   MatDenseGetArray(C1,&u);
 71:   for (i=Istart;i<Iend;i++) {
 72:     if (i<N/2) u[i-Istart] = 1.0;
 73:     if (i==0) u[i+Iend-2*Istart] = -2.0;
 74:     if (i==1) u[i+Iend-2*Istart] = -1.0;
 75:     if (i==2) u[i+Iend-2*Istart] = -1.0;
 76:   }
 77:   MatDenseRestoreArray(C1,&u);
 78:   MatAssemblyBegin(C1,MAT_FINAL_ASSEMBLY);
 79:   MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY);
 80:   MatCreateLRC(NULL,C1,NULL,NULL,&C);
 81:   MatDestroy(&C1);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                 Create the solver and set various options
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 86:   LMECreate(PETSC_COMM_WORLD,&lme);
 87:   LMESetProblemType(lme,LME_SYLVESTER);
 88:   LMEGetProblemType(lme,&ptype);
 89:   PetscPrintf(PETSC_COMM_WORLD," Equation type set to %d\n",ptype);
 90:   LMESetProblemType(lme,LME_LYAPUNOV);
 91:   LMEGetProblemType(lme,&ptype);
 92:   PetscPrintf(PETSC_COMM_WORLD," Equation type changed to %d\n",ptype);
 93:   LMESetCoefficients(lme,A,NULL,NULL,NULL);
 94:   LMESetRHS(lme,C);

 96:   /* test prefix usage */
 97:   if (testprefix) {
 98:     LMESetOptionsPrefix(lme,"check_");
 99:     LMEAppendOptionsPrefix(lme,"myprefix_");
100:     LMEGetOptionsPrefix(lme,&prefix);
101:     PetscPrintf(PETSC_COMM_WORLD," LME prefix is currently: %s\n",prefix);
102:   }

104:   /* test some interface functions */
105:   LMEGetCoefficients(lme,&B,NULL,NULL,NULL);
106:   if (viewmatrices) MatView(B,PETSC_VIEWER_STDOUT_WORLD);
107:   LMEGetRHS(lme,&D);
108:   if (viewmatrices) MatView(D,PETSC_VIEWER_STDOUT_WORLD);
109:   LMESetTolerances(lme,PETSC_DEFAULT,100);
110:   LMESetDimensions(lme,21);
111:   LMESetErrorIfNotConverged(lme,PETSC_TRUE);
112:   /* test monitors */
113:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
114:   LMEMonitorSet(lme,(PetscErrorCode (*)(LME,PetscInt,PetscReal,void*))LMEMonitorDefault,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
115:   /* LMEMonitorCancel(lme); */
116:   LMESetFromOptions(lme);

118:   LMEGetType(lme,&type);
119:   PetscPrintf(PETSC_COMM_WORLD," Solver being used: %s\n",type);

121:   /* query properties and print them */
122:   LMEGetTolerances(lme,&tol,&maxit);
123:   PetscPrintf(PETSC_COMM_WORLD," Tolerance: %g, max iterations: %" PetscInt_FMT "\n",(double)tol,maxit);
124:   LMEGetDimensions(lme,&ncv);
125:   PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv);
126:   LMEGetErrorIfNotConverged(lme,&flg);
127:   if (flg) PetscPrintf(PETSC_COMM_WORLD," Erroring out if convergence fails\n");

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130:                 Solve the matrix equation and compute residual error
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

133:   LMESolve(lme);
134:   LMEGetErrorEstimate(lme,&errest);
135:   PetscPrintf(PETSC_COMM_WORLD," Error estimate reported by the solver: %.4g\n",(double)errest);
136:   LMEComputeError(lme,&error);
137:   PetscPrintf(PETSC_COMM_WORLD," Computed residual norm: %.4g\n\n",(double)error);

139:   /*
140:      Free work space
141:   */
142:   LMEDestroy(&lme);
143:   MatDestroy(&A);
144:   MatDestroy(&C);
145:   SlepcFinalize();
146:   return 0;
147: }

149: /*TEST

151:    test:
152:       suffix: 1
153:       args: -lme_monitor_cancel -lme_converged_reason -lme_view -view_matrices -log_exclude lme,bv
154:       requires: double
155:       filter: sed -e "s/4.0[0-9]*e-10/4.03e-10/"

157:    test:
158:       suffix: 2
159:       args: -test_prefix -check_myprefix_lme_monitor
160:       requires: double
161:       filter: sed -e "s/estimate [0-9]\.[0-9]*e[+-]\([0-9]*\)/estimate (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/"

163:    test:
164:       suffix: 3
165:       args: -lme_monitor_cancel -info -lme_monitor draw::draw_lg -draw_virtual
166:       requires: x double
167:       filter: sed -e "s/equation = [0-9]\.[0-9]*e[+-]\([0-9]*\)/equation = (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/" | grep -v Comm | grep -v machine | grep -v PetscGetHostName | grep -v OpenMP | grep -v Colormap | grep -v "Rank of the Cholesky factor" | grep -v "potrf failed" | grep -v "querying" | grep -v FPTrap | grep -v Device

169: TEST*/