Actual source code: test1.c

slepc-3.19.0 2023-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 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:   PetscFunctionBeginUser;
 32:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

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

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

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

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

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

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

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

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

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

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

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

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

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

150: /*TEST

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

158:    test:
159:       suffix: 2
160:       args: -test_prefix -check_myprefix_lme_monitor
161:       requires: double
162:       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/"

164:    test:
165:       suffix: 3
166:       args: -lme_monitor_cancel -info -lme_monitor draw::draw_lg -draw_virtual
167:       requires: x double
168:       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

170: TEST*/