]> SALOME platform Git repositories - tools/solverlab.git/commitdiff
Salome HOME
Added computation of singular vectors in ProblemCoreFlows
authormichael <michael@localhost.localdomain>
Sat, 1 May 2021 16:44:29 +0000 (18:44 +0200)
committermichael <michael@localhost.localdomain>
Sat, 1 May 2021 16:44:29 +0000 (18:44 +0200)
CoreFlows/Models/inc/ProblemCoreFlows.hxx
CoreFlows/Models/src/ProblemCoreFlows.cxx

index af168699870e23e768597927ec81230555fe29a8..be45af38e8d4ff5c0f99f2702b731d8d78a604d6 100755 (executable)
@@ -116,7 +116,7 @@ public :
 
        /** \fn solveTimeStep
         * \brief calcule les valeurs inconnues au pas de temps +1 .
-        *  \details c'est une fonction virtuelle ,
+        *  \details c'est une fonction virtuelle
         *  @param  void
         *  \return Renvoie false en cas de problème durant le calcul (valeurs non physiques..)
         *  */
@@ -483,6 +483,21 @@ public :
                _path=resultsPath;
        }
 
+       /** \fn getTimeScheme
+        * \brief returns the  time scheme name
+        * \param [in] void
+        * \param [out] enum TimeScheme (explicit or implicit)
+        *  */
+       TimeScheme getTimeScheme();
+
+       /** \fn setNumericalScheme
+        * \brief sets the numerical method ( explicit vs implicit )
+        * \details
+        * \param [in] TimeScheme
+        * \param [out] void
+        *  */
+       void setTimeScheme( TimeScheme method);
+
        //Couplages Thermohydraulique-thermique-neutronique *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
 
        /** \fn setHeatPowerField
@@ -582,11 +597,13 @@ public :
                _system = system;
        };
 
-    //Spectrum analysis
+    //Spectral analysis
     double getConditionNumber(bool isSingular=false, double tol=1e-6) const;
     std::vector< double > getEigenvalues (int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
     std::vector< Vector > getEigenvectors(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
     Field getEigenvectorsField(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
+    std::vector< double > getSingularValues( int nsv, SVDWhich which=SVD_SMALLEST, double tol=1e-6) const;
+    std::vector< Vector > getSingularVectors(int nsv, SVDWhich which=SVD_SMALLEST, double tol=1e-6) const;
 
        //  some supplementary functions
 
@@ -597,7 +614,7 @@ public :
         * @param name, string, name or description of the matrix
         * @return displays the matrix on the terminal
         *  */
-       void displayMatrix(double *matrix, int size, string name);
+       static void displayMatrix(double *matrix, int size, string name="Vector coefficients :");
 
        /** \fn displayMatrix
         * \brief displays a vector of size "size" for profiling
@@ -606,23 +623,7 @@ public :
         * @param name, string, name or description of the vector
         * @return displays the vector on the terminal
         *  */
-       void displayVector(double *vector, int size, string name);
-
-       /** \fn getTimeScheme
-        * \brief returns the  time scheme name
-        * \param [in] void
-        * \param [out] enum TimeScheme (explicit or implicit)
-        *  */
-       TimeScheme getTimeScheme();
-
-       /** \fn setNumericalScheme
-        * \brief sets the numerical method ( explicit vs implicit )
-        * \details
-        * \param [in] TimeScheme
-        * \param [out] void
-        *  */
-       void setTimeScheme( TimeScheme method);
-
+       static void displayVector(double *vector, int size, string name="Matrix coefficients :");
 
 protected :
 
index 11f7e25767c84d42cf951f6f03ad000c732c192b..7d69103ffa60b3aa92db4333326a271a7bae81cd 100755 (executable)
@@ -622,6 +622,19 @@ ProblemCoreFlows::getEigenvectorsField(int nev, EPSWhich which, double tol) cons
   return my_eigenfield;
 }
 
+std::vector< double > 
+ProblemCoreFlows::getSingularValues( int nsv, SVDWhich which, double tol) const
+{
+  SparseMatrixPetsc A = SparseMatrixPetsc(_A);
+  return A.getSingularValues( nsv, which, tol);
+}
+std::vector< Vector > 
+ProblemCoreFlows::getSingularVectors(int nsv, SVDWhich which, double tol) const
+{
+  SparseMatrixPetsc A = SparseMatrixPetsc(_A);
+  return A.getSingularVectors( nsv, which, tol);
+}
+
 Field 
 ProblemCoreFlows::getUnknownField() const
 {