|
Teuchos Package Browser (Single Doxygen Collection)
Version of the Day
|
Go to the documentation of this file.
58 #ifdef HAVE_TEUCHOS_COMPLEX
59 #define STYPE std::complex<double>
66 #ifdef HAVE_TEUCHOS_COMPLEX
67 #define SCALARMAX STYPE(10,0)
69 #define SCALARMAX STYPE(10)
72 template<
typename TYPE>
82 template<
typename TYPE>
94 std::complex<T>
GetRandom( std::complex<T>, std::complex<T> );
107 int main(
int argc,
char* argv[])
115 if (argc>1)
if (argv[1][0]==
'-' && argv[1][1]==
'v') verbose =
true;
120 int numberFailedTests = 0, failedComparison = 0;
122 std::string testName =
"", testType =
"";
124 #ifdef HAVE_TEUCHOS_COMPLEX
125 testType =
"COMPLEX";
130 if (verbose) std::cout<<std::endl<<
"********** CHECKING TEUCHOS SERIAL SPD DENSE SOLVER - " << testType <<
"-VALUED **********"<<std::endl<<std::endl;
138 for (
int j=0; j<
n; ++j) {
139 for (
int i=0; i<=j; ++i) {
140 (*A1full)(i,j) = (*A1)(i,j);
143 for (
int j=0; j<
n; ++j) {
144 for (
int i=j+1; i<
n; ++i) {
151 testName =
"Generating right-hand side vector using A*x, where x is a random vector:";
152 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
165 returnCode = solver1.
factor();
166 testName =
"Simple solve: factor() random A:";
167 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
170 returnCode = solver1.
solve();
171 testName =
"Simple solve: solve() random A (NO_TRANS):";
173 if(verbose && failedComparison>0) std::cout <<
"COMPARISON FAILED : ";
174 numberFailedTests += failedComparison;
175 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
178 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
186 for (
int j=0; j<
n; ++j) {
187 for (
int i=0; i<=j; ++i) {
188 (*A2full)(i,j) = (*A2)(i,j);
191 for (
int j=0; j<
n; ++j) {
192 for (
int i=j+1; i<
n; ++i) {
210 returnCode = solver2.
factor();
211 testName =
"Solve with iterative refinement: factor() random A:";
212 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
215 returnCode = solver2.
solve();
216 testName =
"Solve with iterative refinement: solve() random A (NO_TRANS):";
218 if(verbose && failedComparison>0) std::cout <<
"COMPARISON FAILED : ";
219 numberFailedTests += failedComparison;
220 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
231 for (
int j=0; j<
n; ++j) {
232 for (
int i=0; i<=j; ++i) {
233 (*A3full)(i,j) = (*A3)(i,j);
236 for (
int j=0; j<
n; ++j) {
237 for (
int i=j+1; i<
n; ++i) {
258 returnCode = solver3.
factor();
259 testName =
"Solve with matrix equilibration: factor() random A:";
260 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
263 returnCode = solver3.
solve();
264 testName =
"Solve with matrix equilibration: solve() random A (NO_TRANS):";
266 if(verbose && failedComparison>0) std::cout <<
"COMPARISON FAILED : ";
267 numberFailedTests += failedComparison;
268 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
273 returnCode = solver3.
solve();
274 testName =
"Solve with matrix equilibration: solve() without factor() random A (NO_TRANS):";
276 if(verbose && failedComparison>0) std::cout <<
"COMPARISON FAILED : ";
277 numberFailedTests += failedComparison;
278 numberFailedTests +=
ReturnCodeCheck(testName, returnCode, 0, verbose);
283 if(numberFailedTests > 0)
286 std::cout <<
"Number of failed tests: " << numberFailedTests << std::endl;
287 std::cout <<
"End Result: TEST FAILED" << std::endl;
291 if(numberFailedTests == 0)
292 std::cout <<
"End Result: TEST PASSED" << std::endl;
298 template<
typename TYPE>
299 int PrintTestResults(std::string testName, TYPE calculatedResult, TYPE expectedResult,
bool verbose)
302 if(calculatedResult == expectedResult)
304 if(verbose) std::cout << testName <<
" successful." << std::endl;
309 if(verbose) std::cout << testName <<
" unsuccessful." << std::endl;
315 int ReturnCodeCheck(std::string testName,
int returnCode,
int expectedResult,
bool verbose)
318 if(expectedResult == 0)
322 if(verbose) std::cout << testName <<
" test successful." << std::endl;
327 if(verbose) std::cout << testName <<
" test unsuccessful. Return code was " << returnCode <<
"." << std::endl;
335 if(verbose) std::cout << testName <<
" test successful -- failed as expected." << std::endl;
340 if(verbose) std::cout << testName <<
" test unsuccessful -- did not fail as expected. Return code was " << returnCode <<
"." << std::endl;
347 template<
typename TYPE>
354 std::complex<T>
GetRandom( std::complex<T> Low, std::complex<T> High)
356 T lowMag = Low.real();
357 T highMag = High.real();
360 return std::complex<T>( real, imag );
394 for (
int i=0; i<
n; ++i) {
404 for (
int j=0; j<
n; ++j) {
405 for (
int i=0; i<=j; ++i) {
406 (*mat)(i,j) = (*tmp2)(i,j);
418 for (
int i=0; i<m; i++)
419 for (
int j=0; j<
n; j++)
430 for (
int i=0; i<
n; i++)
450 MagnitudeType temp = norm_diff;
454 if (temp > Tolerance)
int factor()
Computes the in-place QR factorization of the matrix using the LAPACK routine _GETRF or the Eigen cla...
Templated class for constructing and using Hermitian positive definite dense matrices.
Reference-counted pointer class and non-member templated function implementations.
Templated class for solving dense linear problems.
int CompareVectors(const SerialDenseVector< OTYPE, STYPE > &Vector1, const SerialDenseVector< OTYPE, STYPE > &Vector2, ScalarTraits< STYPE >::magnitudeType Tolerance)
int main(int argc, char *argv[])
int ReturnCodeCheck(std::string, int, int, bool)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
void factorWithEquilibration(bool flag)
Causes equilibration to be called just before the matrix factorization as part of the call to factor.
int formQ()
Explicitly forms the unitary matrix Q.
int PrintTestResults(std::string, TYPE, TYPE, bool)
Teuchos::RCP< DVector > GetRandomVector(int n)
This class creates and provides basic support for symmetric, positive-definite dense matrices of temp...
int solve()
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()....
Teuchos::RCP< DMatrix > GetRandomMatrix(int m, int n)
A class for solving dense linear problems.
ScalarTraits< ScalarType >::magnitudeType normFrobenius() const
Returns the Frobenius-norm of the matrix.
Smart reference counting pointer class for automatic garbage collection.
void setUpper()
Specify that the upper triangle of the this matrix should be used.
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)
Sets the pointers for coefficient matrix.
int factor()
Computes the in-place Cholesky factorization of the matrix using the LAPACK routine DPOTRF.
This structure defines some basic traits for a scalar field type.
A class for constructing and using Hermitian positive definite dense matrices.
SerialDenseVector< OTYPE, STYPE > DVector
Templated serial dense matrix class.
TYPE GetRandom(TYPE, TYPE)
SerialDenseMatrix< OTYPE, STYPE > DMatrix
std::string Teuchos_Version()
This class creates and provides basic support for dense vectors of templated type as a specialization...
Templated serial, dense, symmetric matrix class.
This class creates and provides basic support for dense rectangular matrix of templated type.
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
Multiply A * B and add them to this; this = beta * this + alpha*A*B.
SerialSymDenseMatrix< OTYPE, STYPE > SDMatrix
Defines basic traits for the scalar field type.
Teuchos::RCP< SDMatrix > GetRandomSpdMatrix(int n)
int setMatrix(const RCP< SerialSymDenseMatrix< OrdinalType, ScalarType > > &A_in)
Sets the pointers for coefficient matrix.
Templated serial dense vector class.
Non-member helper functions on the templated serial, dense matrix/vector classes.
void solveToRefinedSolution(bool flag)
Causes all solves to compute solution to best ability using iterative refinement.
int setVectors(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &X, const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &B)
Sets the pointers for left and right hand side vector(s).
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getQ() const
Returns pointer to Q (assuming factorization has been performed).