Intrepid2
Intrepid2_Utils_ExtDataDef.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
15 #ifndef __INTREPID2_UTILS_EXTDATA_DEF_HPP__
16 #define __INTREPID2_UTILS_EXTDATA_DEF_HPP__
17 
18 namespace Intrepid2 {
19 
20  /***************************************************************************************************
21  * *
22  * Utility functions for handling external data in tests *
23  * *
24  ***************************************************************************************************/
25 
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");
35 
36  Teuchos::RCP<std::ostream> outStream;
37  Teuchos::oblackholestream outNothing;
38  if (iprint > 0)
39  outStream = Teuchos::rcp(&std::cout, false);
40  else
41  outStream = Teuchos::rcp(&outNothing, false);
42 
43  // Save the format state of the original std::cout.
44  Teuchos::oblackholestream oldFormatState;
45  oldFormatState.copyfmt(std::cout);
46 
47  std::string line;
48  ValueType testentry;
49  ValueType abstol;
50  ValueType absdiff;
51  ordinal_type i = 0, j = 0;
52  ordinal_type err = 0;
53 
54  while (! inputFile.eof() && i < static_cast<ordinal_type>(testMat.extent(0)) ) {
55  std::getline(inputFile,line);
56  std::istringstream linestream(line);
57  std::string chunk;
58  j = 0;
59  while( linestream >> chunk ) {
60  ordinal_type num1;
61  ordinal_type num2;
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);
66  chunkstream >> num1;
67  chunkstream >> num2;
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) {
72  ++err;
73  *outStream << "FAILURE --> ";
74  }
75  *outStream << "entry[" << i << "," << j << "]:" << " "
76  << testMat(i, j) << " " << num1 << "/" << num2 << " "
77  << absdiff << " " << "<?" << " " << abstol << "\n";
78  }
79  else {
80  std::istringstream chunkstream(chunk);
81  if (analyticDataType == INTREPID2_UTILS_FRACTION) {
82  chunkstream >> num1;
83  testentry = (ValueType)(num1);
84  }
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) {
90  ++err;
91  *outStream << "FAILURE --> ";
92  }
93  *outStream << "entry[" << i << "," << j << "]:" << " "
94  << testMat(i, j) << " " << testentry << " "
95  << absdiff << " " << "<?" << " " << abstol << "\n";
96  }
97  ++j;
98  }
99  ++i;
100  }
101 
102  // reset format state of std::cout
103  std::cout.copyfmt(oldFormatState);
104 
105  return err;
106  } // end compareToAnalytic
107 
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");
115 
116  Teuchos::oblackholestream oldFormatState;
117  oldFormatState.copyfmt(std::cout);
118 
119  std::string line;
120  ValueType testentry;
121  ordinal_type i = 0, j = 0;
122 
123  while (! inputFile.eof() && i < static_cast<ordinal_type>(testMat.extent(0)) ) {
124  std::getline (inputFile,line);
125  std::istringstream linestream(line);
126  std::string chunk;
127  j = 0;
128  while( linestream >> chunk ) {
129  ordinal_type num1;
130  ordinal_type num2;
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);
135  chunkstream >> num1;
136  chunkstream >> num2;
137  testentry = (ValueType)(num1)/(ValueType)(num2);
138  testMat(i, j) = testentry;
139  }
140  else {
141  std::istringstream chunkstream(chunk);
142  if (analyticDataType == INTREPID2_UTILS_FRACTION) {
143  chunkstream >> num1;
144  testentry = (ValueType)(num1);
145  }
146  else if (analyticDataType == INTREPID2_UTILS_SCALAR)
147  chunkstream >> testentry;
148  testMat(i, j) = testentry;
149  }
150  ++j;
151  }
152  ++i;
153  }
154 
155  // reset format state of std::cout
156  std::cout.copyfmt(oldFormatState);
157  } // end getAnalytic
158 
159 } // end namespace Intrepid2
160 
161 #endif