]> SALOME platform Git repositories - tools/solverlab.git/commitdiff
Salome HOME
Updated SLEPc syntax.\n Corrected division by zero problem.
authormichael <michael@localhost.localdomain>
Sat, 30 Oct 2021 20:33:43 +0000 (22:33 +0200)
committermichael <michael@localhost.localdomain>
Sat, 30 Oct 2021 20:33:43 +0000 (22:33 +0200)
CDMATH/linearsolver/src/SparseMatrixPetsc.cxx

index d9cd0f934f32e0b9dff1de06c58599d2331ff84c..91412cdf858b17cf3edddbecdb8824ba7b42cfb6 100755 (executable)
@@ -625,13 +625,23 @@ SparseMatrixPetsc::computeSpectrum(int nev, double ** valP, double ***vecP, EPSW
       /*
          Compute the relative error associated to each eigenpair
       */
-      EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);
-
-      if (ki!=0.0) {
-        PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)kr,(double)ki,(double)error);
-      } else {
-        PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)kr,(double)error);
-      }
+      if(fabs(kr)>tol || fabs(ki)>tol)
+      {
+             EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);
+       
+             if (ki!=0.0)
+               PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)kr,(double)ki,(double)error);
+             else
+               PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)kr,(double)error);
+          }
+          else
+      {
+             if (ki!=0.0)
+               PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12s\n",(double)kr,(double)ki,"Null eigenvalue");
+             else
+               PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12s\n",(double)kr,"Null eigenvalue");
+          }
+          
       *(*valP + i)=kr;
       VecGetArray(xr,&myVecp);
       *(*vecP+  i)=new double [_numberOfRows];
@@ -691,7 +701,12 @@ SparseMatrixPetsc::computeSVD(int nsv, double ** valS, double ***vecS, SVDWhich
   /*
      Set operators. In this case, it is a standard singular value problem
   */
-  SVDSetOperator(svd,_mat);
+#if (SLEPC_VERSION_MAJOR==3) && (SLEPC_VERSION_MINOR > 14)
+       SVDSetOperators(svd,_mat,NULL);
+#else
+       SVDSetOperator(svd,_mat);
+#endif
+  
   SVDSetWhichSingularTriplets(svd,which);
   SVDSetDimensions(svd,nsv,PETSC_DEFAULT,PETSC_DEFAULT);
   SVDSetTolerances(svd,tol,PETSC_DEFAULT);
@@ -744,13 +759,18 @@ SparseMatrixPetsc::computeSVD(int nsv, double ** valS, double ***vecS, SVDWhich
         Get converged singular values: i-th eigenvalue is stored in valS
       */
       SVDGetSingularTriplet(svd,i,*valS+i, u, v);
-      /*
-         Compute the relative error associated to each singular value
-      */
-      SVDComputeError( svd, i, SVD_ERROR_RELATIVE, &error );
-
-      PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)*(*valS+i),(double)error);
-
+      if(fabs(*(*valS+i))>tol)
+      {
+             /*
+                Compute the relative error associated to each singular value
+             */
+             SVDComputeError( svd, i, SVD_ERROR_RELATIVE, &error );
+       
+             PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)*(*valS+i),(double)error);
+               }
+       else
+          PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12s\n",(double)*(*valS+i),"Null singular value");
+       
       VecGetArray(u,&myVecS);
       *(*vecS + i)=new double [_numberOfRows];
       memcpy(*(*vecS + i),myVecS,_numberOfRows*sizeof(double)) ;