15 #ifndef __INTREPID2_UTILS_EXTDATA_DEF_HPP__
16 #define __INTREPID2_UTILS_EXTDATA_DEF_HPP__
26 template<
typename ValueType,
27 class ...testMatProperties>
28 ordinal_type compareToAnalytic( std::ifstream &inputFile,
29 const Kokkos::DynRankView<ValueType,testMatProperties...> testMat,
30 const ValueType reltol,
31 const ordinal_type iprint,
32 const TypeOfExactData analyticDataType ) {
33 INTREPID2_TEST_FOR_EXCEPTION( testMat.rank() != 2, std::invalid_argument,
34 ">>> ERROR (compareToAnalytic): testMat must have rank 2");
36 Teuchos::RCP<std::ostream> outStream;
37 Teuchos::oblackholestream outNothing;
39 outStream = Teuchos::rcp(&std::cout,
false);
41 outStream = Teuchos::rcp(&outNothing,
false);
44 Teuchos::oblackholestream oldFormatState;
45 oldFormatState.copyfmt(std::cout);
51 ordinal_type i = 0, j = 0;
54 while (! inputFile.eof() && i < static_cast<ordinal_type>(testMat.extent(0)) ) {
55 std::getline(inputFile,line);
56 std::istringstream linestream(line);
59 while( linestream >> chunk ) {
62 std::string::size_type loc = chunk.find(
"/", 0);
63 if( loc != std::string::npos ) {
64 chunk.replace( loc, 1,
" ");
65 std::istringstream chunkstream(chunk);
68 testentry = (ValueType)(num1)/(ValueType)(num2);
69 abstol = ( std::fabs(testentry) < reltol ? reltol : std::fabs(reltol*testentry) );
70 absdiff = std::fabs(testentry - testMat(i, j));
71 if (absdiff > abstol) {
73 *outStream <<
"FAILURE --> ";
75 *outStream <<
"entry[" << i <<
"," << j <<
"]:" <<
" "
76 << testMat(i, j) <<
" " << num1 <<
"/" << num2 <<
" "
77 << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
80 std::istringstream chunkstream(chunk);
81 if (analyticDataType == INTREPID2_UTILS_FRACTION) {
83 testentry = (ValueType)(num1);
85 else if (analyticDataType == INTREPID2_UTILS_SCALAR)
86 chunkstream >> testentry;
87 abstol = ( std::fabs(testentry) < reltol ?reltol : std::fabs(reltol*testentry) );
88 absdiff = std::fabs(testentry - testMat(i, j));
89 if (absdiff > abstol) {
91 *outStream <<
"FAILURE --> ";
93 *outStream <<
"entry[" << i <<
"," << j <<
"]:" <<
" "
94 << testMat(i, j) <<
" " << testentry <<
" "
95 << absdiff <<
" " <<
"<?" <<
" " << abstol <<
"\n";
103 std::cout.copyfmt(oldFormatState);
108 template<
typename ValueType,
109 class ...testMatProperties>
110 void getAnalytic( Kokkos::DynRankView<ValueType,testMatProperties...> testMat,
111 std::ifstream &inputFile,
112 const TypeOfExactData analyticDataType ) {
113 INTREPID2_TEST_FOR_EXCEPTION( testMat.rank() != 2, std::invalid_argument,
114 ">>> ERROR (getToAnalytic): testMat must have rank 2");
116 Teuchos::oblackholestream oldFormatState;
117 oldFormatState.copyfmt(std::cout);
121 ordinal_type i = 0, j = 0;
123 while (! inputFile.eof() && i < static_cast<ordinal_type>(testMat.extent(0)) ) {
124 std::getline (inputFile,line);
125 std::istringstream linestream(line);
128 while( linestream >> chunk ) {
131 std::string::size_type loc = chunk.find(
"/", 0);
132 if( loc != std::string::npos ) {
133 chunk.replace( loc, 1,
" ");
134 std::istringstream chunkstream(chunk);
137 testentry = (ValueType)(num1)/(ValueType)(num2);
138 testMat(i, j) = testentry;
141 std::istringstream chunkstream(chunk);
142 if (analyticDataType == INTREPID2_UTILS_FRACTION) {
144 testentry = (ValueType)(num1);
146 else if (analyticDataType == INTREPID2_UTILS_SCALAR)
147 chunkstream >> testentry;
148 testMat(i, j) = testentry;
156 std::cout.copyfmt(oldFormatState);