00001
00002
00003
00004
00005
00006
00007 #ifndef ALBANY_MPAS_STKMESHSTRUCT_HPP
00008 #define ALBANY_MPAS_STKMESHSTRUCT_HPP
00009
00010 #include "Albany_GenericSTKMeshStruct.hpp"
00011
00012
00013 namespace Albany {
00014
00015 class MpasSTKMeshStruct : public GenericSTKMeshStruct {
00016
00017 public:
00018
00019 MpasSTKMeshStruct(const Teuchos::RCP<Teuchos::ParameterList>& params,
00020 const Teuchos::RCP<const Epetra_Comm>& comm,
00021 const std::vector<int>& indexToTriangleID, const std::vector<int>& verticesOnTria, int nGlobalTriangles);
00022
00023 MpasSTKMeshStruct(const Teuchos::RCP<Teuchos::ParameterList>& params,
00024 const Teuchos::RCP<const Epetra_Comm>& comm,
00025 const std::vector<int>& indexToTriangleID, const std::vector<int>& verticesOnTria, int nGlobalTriangles, int numLayers, int Ordering = 0);
00026
00027 MpasSTKMeshStruct(const Teuchos::RCP<Teuchos::ParameterList>& params,
00028 const Teuchos::RCP<const Epetra_Comm>& comm,
00029 const std::vector<int>& indexToTriangleID, int nGlobalTriangles, int numLayers, int Ordering = 0);
00030
00031
00032 ~MpasSTKMeshStruct();
00033
00034 void setFieldAndBulkData(
00035 const Teuchos::RCP<const Epetra_Comm>& comm,
00036 const Teuchos::RCP<Teuchos::ParameterList>& params,
00037 const unsigned int neq_,
00038 const Albany::AbstractFieldContainer::FieldContainerRequirements& req,
00039 const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00040 const unsigned int worksetSize){};
00041
00043 bool hasRestartSolution() const {return hasRestartSol; }
00044
00045 void setHasRestartSolution(bool hasRestartSolution) {hasRestartSol = hasRestartSolution; }
00046
00047 void setRestartDataTime(double restartT) {restartTime = restartT; }
00048
00050 double restartDataTime() const {return restartTime;}
00051
00052
00053 void constructMesh(
00054 const Teuchos::RCP<const Epetra_Comm>& comm,
00055 const Teuchos::RCP<Teuchos::ParameterList>& params,
00056 const unsigned int neq_,
00057 const Albany::AbstractFieldContainer::FieldContainerRequirements& req,
00058 const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00059 const std::vector<int>& indexToVertexID, const std::vector<double>& verticesCoords, const std::vector<bool>& isVertexBoundary, int nGlobalVertices,
00060 const std::vector<int>& verticesOnTria,
00061 const std::vector<bool>& isBoundaryEdge, const std::vector<int>& trianglesOnEdge, const std::vector<int>& trianglesPositionsOnEdge,
00062 const std::vector<int>& verticesOnEdge, const std::vector<int>& indexToEdgeID, int nGlobalEdges,
00063 const unsigned int worksetSize);
00064
00065 void constructMesh(
00066 const Teuchos::RCP<const Epetra_Comm>& comm,
00067 const Teuchos::RCP<Teuchos::ParameterList>& params,
00068 const unsigned int neq_,
00069 const Albany::AbstractFieldContainer::FieldContainerRequirements& req,
00070 const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00071 const std::vector<int>& indexToVertexID, const std::vector<double>& verticesCoords, const std::vector<bool>& isVertexBoundary, int nGlobalVertices,
00072 const std::vector<int>& verticesOnTria,
00073 const std::vector<bool>& isBoundaryEdge, const std::vector<int>& trianglesOnEdge, const std::vector<int>& trianglesPositionsOnEdge,
00074 const std::vector<int>& verticesOnEdge,
00075 const std::vector<int>& indexToEdgeID, int nGlobalEdges,
00076 const std::vector<int>& indexToTriangleID,
00077 const unsigned int worksetSize,
00078 int numLayers, int Ordering = 0);
00079
00080 void constructMesh(
00081 const Teuchos::RCP<const Epetra_Comm>& comm,
00082 const Teuchos::RCP<Teuchos::ParameterList>& params,
00083 const unsigned int neq_,
00084 const Albany::AbstractFieldContainer::FieldContainerRequirements& req,
00085 const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00086 const std::vector<int>& indexToVertexID, const std::vector<int>& indexToMpasVertexID, const std::vector<double>& verticesCoords, const std::vector<bool>& isVertexBoundary, int nGlobalVertices,
00087 const std::vector<int>& verticesOnTria,
00088 const std::vector<bool>& isBoundaryEdge, const std::vector<int>& trianglesOnEdge, const std::vector<int>& trianglesPositionsOnEdge,
00089 const std::vector<int>& verticesOnEdge,
00090 const std::vector<int>& indexToEdgeID, int nGlobalEdges,
00091 const std::vector<int>& indexToTriangleID,
00092 const unsigned int worksetSize,
00093 int numLayers, int Ordering = 0);
00094
00095
00096 const bool getInterleavedOrdering() const {return this->interleavedOrdering;}
00097
00098 private:
00099
00100 inline void tetrasFromPrismStructured (int const* prismVertexLIds, int const* prismVertexGIds, int tetrasIdsOnPrism[][4]);
00101 inline void setBdFacesOnPrism (const std::vector<std::vector<std::vector<int> > >& prismStruct, const std::vector<int>& prismFaceIds, std::vector<int>& tetraPos, std::vector<int>& facePos);
00102
00103 Teuchos::RCP<const Teuchos::ParameterList>
00104 getValidDiscretizationParameters() const;
00105
00106 Teuchos::RCP<Teuchos::FancyOStream> out;
00107 bool periodic;
00108 int NumEles;
00109 bool hasRestartSol;
00110 double restartTime;
00111 Teuchos::RCP<Epetra_Map> elem_map;
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130 };
00131
00132
00133
00134 inline void MpasSTKMeshStruct::tetrasFromPrismStructured (int const* prismVertexMpasIds, int const* prismVertexGIds, int tetrasIdsOnPrism[][4])
00135 {
00136 int PrismVerticesMap[6][6] = {{0, 1, 2, 3, 4, 5}, {1, 2, 0, 4, 5, 3}, {2, 0, 1, 5, 3, 4}, {3, 5, 4, 0, 2, 1}, {4, 3, 5, 1, 0, 2}, {5, 4, 3, 2, 1, 0}};
00137
00138 int tetraOfPrism[2][3][4] = {{{0, 1, 2, 5}, {0, 1, 5, 4}, {0, 4, 5, 3}}, {{0, 1, 2, 4}, {0, 4, 2, 5}, {0, 4, 5, 3}}};
00139
00140 int minIndex = std::min_element (prismVertexMpasIds, prismVertexMpasIds + 3) - prismVertexMpasIds;
00141
00142 int v1 (prismVertexMpasIds[PrismVerticesMap[minIndex][1]]);
00143 int v2 (prismVertexMpasIds[PrismVerticesMap[minIndex][2]]);
00144
00145 int prismType = v1 > v2;
00146
00147 for (int iTetra = 0; iTetra < 3; iTetra++)
00148 for (int iVertex = 0; iVertex < 4; iVertex++)
00149 {
00150 tetrasIdsOnPrism[iTetra][iVertex] = prismVertexGIds[tetraOfPrism[prismType][iTetra][iVertex]];
00151 }
00152
00153
00154
00155 int reorderedPrismLIds[6];
00156
00157 for (int ii = 0; ii < 6; ii++)
00158 {
00159 reorderedPrismLIds[ii] = prismVertexGIds[PrismVerticesMap[minIndex][ii]];
00160 }
00161
00162 for (int iTetra = 0; iTetra < 3; iTetra++)
00163 for (int iVertex = 0; iVertex < 4; iVertex++)
00164 {
00165 tetrasIdsOnPrism[iTetra][iVertex] = reorderedPrismLIds[tetraOfPrism[prismType][iTetra][iVertex]];
00166 }
00167 }
00168
00169
00170
00171
00172 void MpasSTKMeshStruct::setBdFacesOnPrism (const std::vector<std::vector<std::vector<int> > >& prismStruct, const std::vector<int>& prismFaceIds, std::vector<int>& tetraPos, std::vector<int>& facePos)
00173 {
00174 int numTriaFaces = prismFaceIds.size() - 2;
00175 tetraPos.assign(numTriaFaces,-1);
00176 facePos.assign(numTriaFaces,-1);
00177
00178
00179 for (int iTetra (0), k (0); (iTetra < 3 && k < numTriaFaces); iTetra++)
00180 {
00181 bool found;
00182 for (int jFaceLocalId = 0; jFaceLocalId < 4; jFaceLocalId++ )
00183 {
00184 found = true;
00185 for (int ip (0); ip < 3 && found; ip++)
00186 {
00187 int localId = prismStruct[iTetra][jFaceLocalId][ip];
00188 int j = 0;
00189 found = false;
00190 while ( (j < prismFaceIds.size()) && !found )
00191 {
00192 found = (localId == prismFaceIds[j]);
00193 j++;
00194 }
00195 }
00196 if (found)
00197 {
00198 tetraPos[k] = iTetra;
00199 facePos[k] = jFaceLocalId;
00200 k += found;
00201 break;
00202 }
00203 }
00204 }
00205 }
00206
00207 }
00208
00209 #endif