00001
00002
00003
00004
00005
00006
00007 #ifndef ALBPUMI_QPDATA_HPP
00008 #define ALBPUMI_QPDATA_HPP
00009
00010
00011 #include "Teuchos_RCP.hpp"
00012 #include "Teuchos_ParameterList.hpp"
00013
00014 #include "PHAL_Dimension.hpp"
00015 #include "Albany_StateInfoStruct.hpp"
00016
00017 namespace AlbPUMI {
00018
00019
00020 template<typename DataType, unsigned Dim>
00021 struct QPData_Traits { };
00022
00023 template<typename DataType, unsigned Dim, class traits = QPData_Traits<DataType, Dim> >
00024 class QPData {
00025
00026 public:
00027
00028 QPData(const std::string& name, const std::vector<int>& dim, const bool output = false);
00029 ~QPData(){}
00030
00032 typedef traits traits_type;
00033
00035 typedef typename traits_type::field_type field_type;
00036
00037 void reAllocateBuffer(const std::size_t nelems);
00038 Albany::MDArray getMDA(const std::size_t nElemsInBucket);
00039
00040 const std::string name;
00041 const bool output;
00042 std::vector<DataType> buffer;
00043 std::vector<int> dims;
00044 int nfield_dofs;
00045 std::size_t beginning_index;
00046
00047 };
00048
00049
00050
00051
00052 template <typename T>
00053 struct QPData_Traits<T, 1> {
00054
00055 enum { size = 1 };
00056 typedef shards::Array<T, shards::NaturalOrder, Cell> field_type ;
00057 static Albany::MDArray buildArray(T *buf, unsigned nelems, std::vector<int>& dims){
00058
00059 return field_type(buf, nelems);
00060
00061 }
00062
00063 };
00064
00065
00066 template <typename T>
00067 struct QPData_Traits<T, 2> {
00068
00069 enum { size = 2 };
00070 typedef shards::Array<T, shards::NaturalOrder, Cell, QuadPoint> field_type ;
00071 static Albany::MDArray buildArray(T *buf, unsigned nelems, std::vector<int>& dims){
00072
00073 return field_type(buf, nelems, dims[1]);
00074
00075 }
00076
00077 };
00078
00079
00080 template <typename T>
00081 struct QPData_Traits<T, 3> {
00082
00083 enum { size = 3 };
00084 typedef shards::Array<T, shards::NaturalOrder, Cell, QuadPoint, Dim> field_type ;
00085 static Albany::MDArray buildArray(T *buf, unsigned nelems, std::vector<int>& dims){
00086
00087 return field_type(buf, nelems, dims[1], dims[2]);
00088
00089 }
00090
00091 };
00092
00093
00094 template <typename T>
00095 struct QPData_Traits<T, 4> {
00096
00097 enum { size = 4 };
00098 typedef shards::Array<T, shards::NaturalOrder, Cell, QuadPoint, Dim, Dim> field_type ;
00099 static Albany::MDArray buildArray(T *buf, unsigned nelems, std::vector<int>& dims){
00100
00101 return field_type(buf, nelems, dims[1], dims[2], dims[3]);
00102
00103 }
00104
00105 };
00106
00107
00108 }
00109
00110
00111 #define QPDATA_INSTANTIATE_TEMPLATE_CLASS_1D(type, name) \
00112 template class name<type, 1>;
00113 #define QPDATA_INSTANTIATE_TEMPLATE_CLASS_2D(type, name) \
00114 template class name<type, 2>;
00115 #define QPDATA_INSTANTIATE_TEMPLATE_CLASS_3D(type, name) \
00116 template class name<type, 3>;
00117 #define QPDATA_INSTANTIATE_TEMPLATE_CLASS_4D(type, name) \
00118 template class name<type, 4>;
00119
00120 #define QPDATA_INSTANTIATE_TEMPLATE_CLASS(name) \
00121 QPDATA_INSTANTIATE_TEMPLATE_CLASS_1D(double, name) \
00122 QPDATA_INSTANTIATE_TEMPLATE_CLASS_2D(double, name) \
00123 QPDATA_INSTANTIATE_TEMPLATE_CLASS_3D(double, name) \
00124 QPDATA_INSTANTIATE_TEMPLATE_CLASS_4D(double, name)
00125
00126 #endif