48 #ifndef __INTREPID2_UTILS_EXTDATA_DEF_HPP__
49 #define __INTREPID2_UTILS_EXTDATA_DEF_HPP__
59 template<
typename ValueType,
60 class ...testMatProperties>
61 ordinal_type compareToAnalytic( std::ifstream &inputFile,
62 const Kokkos::DynRankView<ValueType,testMatProperties...> testMat,
63 const ValueType reltol,
64 const ordinal_type iprint,
65 const TypeOfExactData analyticDataType ) {
66 INTREPID2_TEST_FOR_EXCEPTION( testMat.rank() != 2, std::invalid_argument,
67 ">>> ERROR (compareToAnalytic): testMat must have rank 2");
69 Teuchos::RCP<std::ostream> outStream;
70 Teuchos::oblackholestream outNothing;
72 outStream = Teuchos::rcp(&std::cout,
false);
74 outStream = Teuchos::rcp(&outNothing,
false);
77 Teuchos::oblackholestream oldFormatState;
78 oldFormatState.copyfmt(std::cout);
84 ordinal_type i = 0, j = 0;
87 while (! inputFile.eof() && i < static_cast<ordinal_type>(testMat.extent(0)) ) {
88 std::getline(inputFile,line);
89 std::istringstream linestream(line);
92 while( linestream >> chunk ) {
95 std::string::size_type loc = chunk.find(
"/", 0);
96 if( loc != std::string::npos ) {
97 chunk.replace( loc, 1,
" ");
98 std::istringstream chunkstream(chunk);
101 testentry = (ValueType)(num1)/(ValueType)(num2);
102 abstol = ( std::fabs(testentry) < reltol ? reltol : std::fabs(reltol*testentry) );
103 absdiff = std::fabs(testentry - testMat(i, j));
104 if (absdiff > abstol) {
106 *outStream <<
"FAILURE --> ";
108 *outStream <<
"entry[" << i <<
"," << j <<
"]:" <<
" "
109 << testMat(i, j) <<
" " << num1 <<
"/" << num2 <<
" "
110 << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
113 std::istringstream chunkstream(chunk);
114 if (analyticDataType == INTREPID2_UTILS_FRACTION) {
116 testentry = (ValueType)(num1);
118 else if (analyticDataType == INTREPID2_UTILS_SCALAR)
119 chunkstream >> testentry;
120 abstol = ( std::fabs(testentry) < reltol ?reltol : std::fabs(reltol*testentry) );
121 absdiff = std::fabs(testentry - testMat(i, j));
122 if (absdiff > abstol) {
124 *outStream <<
"FAILURE --> ";
126 *outStream <<
"entry[" << i <<
"," << j <<
"]:" <<
" "
127 << testMat(i, j) <<
" " << testentry <<
" "
128 << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
136 std::cout.copyfmt(oldFormatState);
141 template<
typename ValueType,
142 class ...testMatProperties>
143 void getAnalytic( Kokkos::DynRankView<ValueType,testMatProperties...> testMat,
144 std::ifstream &inputFile,
145 const TypeOfExactData analyticDataType ) {
146 INTREPID2_TEST_FOR_EXCEPTION( testMat.rank() != 2, std::invalid_argument,
147 ">>> ERROR (getToAnalytic): testMat must have rank 2");
149 Teuchos::oblackholestream oldFormatState;
150 oldFormatState.copyfmt(std::cout);
154 ordinal_type i = 0, j = 0;
156 while (! inputFile.eof() && i < static_cast<ordinal_type>(testMat.extent(0)) ) {
157 std::getline (inputFile,line);
158 std::istringstream linestream(line);
161 while( linestream >> chunk ) {
164 std::string::size_type loc = chunk.find(
"/", 0);
165 if( loc != std::string::npos ) {
166 chunk.replace( loc, 1,
" ");
167 std::istringstream chunkstream(chunk);
170 testentry = (ValueType)(num1)/(ValueType)(num2);
171 testMat(i, j) = testentry;
174 std::istringstream chunkstream(chunk);
175 if (analyticDataType == INTREPID2_UTILS_FRACTION) {
177 testentry = (ValueType)(num1);
179 else if (analyticDataType == INTREPID2_UTILS_SCALAR)
180 chunkstream >> testentry;
181 testMat(i, j) = testentry;
189 std::cout.copyfmt(oldFormatState);