Panzer  Version of the Day
Panzer_STK_Utilities.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #include "PanzerAdaptersSTK_config.hpp"
44 
45 #include "Panzer_STK_Utilities.hpp"
47 
48 #include "Kokkos_DynRankView.hpp"
49 
50 #include <stk_mesh/base/FieldBase.hpp>
51 
52 namespace panzer_stk {
53 
54 template <typename GlobalOrdinal>
55 static void gather_in_block(const std::string & blockId, const panzer::UniqueGlobalIndexer<int,GlobalOrdinal> & dofMngr,
56  const Epetra_Vector & x,const std::vector<std::size_t> & localCellIds,
57  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc);
58 
59 #ifdef PANZER_HAVE_FEI
60 void scatter_to_vector(const std::string & blockId, const panzer::DOFManagerFEI<int,int> & dofMngr,
61  const std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc,
62  const std::vector<std::size_t> & localCellIds,
63  Epetra_Vector & x);
64 #endif
65 
66 static void build_local_ids(const panzer_stk::STK_Interface & mesh,
67  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > & localIds);
68 
69 void write_cell_data(panzer_stk::STK_Interface & mesh,const std::vector<double> & data,const std::string & fieldName)
70 {
71  std::vector<std::string> blocks;
72  mesh.getElementBlockNames(blocks);
73 
74  // loop over element blocks
75  for(std::size_t eb=0;eb<blocks.size();eb++) {
76  const std::string & blockId = blocks[eb];
78 
79  std::vector<stk::mesh::Entity> elements;
80  mesh.getMyElements(blockId,elements);
81 
82  // loop over elements in this block
83  for(std::size_t el=0;el<elements.size();el++) {
84  std::size_t localId = mesh.elementLocalId(elements[el]);
85  double * solnData = stk::mesh::field_data(*field,elements[el]);
86  TEUCHOS_ASSERT(solnData!=0); // sanity check
87  solnData[0] = data[localId];
88  }
89  }
90 }
91 
92 template <typename GlobalOrdinal>
93 void write_solution_data(const panzer::UniqueGlobalIndexer<int,GlobalOrdinal> & dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_MultiVector & x,const std::string & prefix,const std::string & postfix)
94 {
95  write_solution_data(dofMngr,mesh,*x(0),prefix,postfix);
96 }
97 
98 template
99 void write_solution_data<int>(const panzer::UniqueGlobalIndexer<int,int> & dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_MultiVector & x,const std::string & prefix,const std::string & postfix);
100 #ifdef PANZER_HAVE_LONG_LONG_INT
101 template
102 void write_solution_data<panzer::Ordinal64>(const panzer::UniqueGlobalIndexer<int,panzer::Ordinal64> & dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_MultiVector & x,const std::string & prefix,const std::string & postfix);
103 #endif
104 
105 template <typename GlobalOrdinal>
106 void write_solution_data(const panzer::UniqueGlobalIndexer<int,GlobalOrdinal> & dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_Vector & x,const std::string & prefix,const std::string & postfix)
107 {
108  typedef Kokkos::DynRankView<double,PHX::Device> FieldContainer;
109 
110  // get local IDs
111  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > localIds;
112  build_local_ids(mesh,localIds);
113 
114  // loop over all element blocks
115  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > >::const_iterator itr;
116  for(itr=localIds.begin();itr!=localIds.end();++itr) {
117  std::string blockId = itr->first;
118  const std::vector<std::size_t> & localCellIds = *(itr->second);
119 
120  std::map<std::string,FieldContainer> data;
121 
122  // get all solution data for this block
123  gather_in_block(blockId,dofMngr,x,localCellIds,data);
124 
125  // write out to stk mesh
126  std::map<std::string,FieldContainer>::iterator dataItr;
127  for(dataItr=data.begin();dataItr!=data.end();++dataItr)
128  mesh.setSolutionFieldData(prefix+dataItr->first+postfix,blockId,localCellIds,dataItr->second);
129  }
130 }
131 
132 template
133 void write_solution_data<int>(const panzer::UniqueGlobalIndexer<int,int> & dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_Vector & x,const std::string & prefix,const std::string & postfix);
134 #ifdef PANZER_HAVE_LONG_LONG_INT
135 template
136 void write_solution_data<panzer::Ordinal64>(const panzer::UniqueGlobalIndexer<int,panzer::Ordinal64> & dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_Vector & x,const std::string & prefix,const std::string & postfix);
137 #endif
138 
139 #ifdef PANZER_HAVE_FEI
140 void read_solution_data(const panzer::DOFManagerFEI<int,int> & dofMngr,const panzer_stk::STK_Interface & mesh,Epetra_MultiVector & x)
141 {
142  read_solution_data(dofMngr,mesh,*x(0));
143 }
144 
145 void read_solution_data(const panzer::DOFManagerFEI<int,int> & dofMngr,const panzer_stk::STK_Interface & mesh,Epetra_Vector & x)
146 {
147  typedef Kokkos::DynRankView<double,PHX::Device> FieldContainer;
148 
149  // get local IDs
150  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > localIds;
151  build_local_ids(mesh,localIds);
152 
153  // loop over all element blocks
154  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > >::const_iterator itr;
155  for(itr=localIds.begin();itr!=localIds.end();++itr) {
156  std::string blockId = itr->first;
157  const std::vector<std::size_t> & localCellIds = *(itr->second);
158 
159  std::map<std::string,FieldContainer> data;
160  const std::set<int> & fieldNums = dofMngr.getFields(blockId);
161 
162  // write out to stk mesh
163  std::set<int>::const_iterator fieldItr;
164  for(fieldItr=fieldNums.begin();fieldItr!=fieldNums.end();++fieldItr) {
165  std::string fieldStr = dofMngr.getFieldString(*fieldItr);
166  mesh.getSolutionFieldData(fieldStr,blockId,localCellIds,data[fieldStr]);
167  }
168 
169  // get all solution data for this block
170  scatter_to_vector(blockId,dofMngr,data,localCellIds,x);
171  }
172 }
173 #endif
174 
175 template <typename GlobalOrdinal>
176 void gather_in_block(const std::string & blockId, const panzer::UniqueGlobalIndexer<int,GlobalOrdinal> & dofMngr,
177  const Epetra_Vector & x,const std::vector<std::size_t> & localCellIds,
178  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc)
179 {
180  const std::vector<int> & fieldNums = dofMngr.getBlockFieldNumbers(blockId);
181 
182  for(std::size_t fieldIndex=0;fieldIndex<fieldNums.size();fieldIndex++) {
183  int fieldNum = fieldNums[fieldIndex];
184  std::string fieldStr = dofMngr.getFieldString(fieldNum);
185 
186  // grab the field
187  const std::vector<int> & elmtOffset = dofMngr.getGIDFieldOffsets(blockId,fieldNum);
188  fc[fieldStr] = Kokkos::DynRankView<double,PHX::Device>("fc",localCellIds.size(),elmtOffset.size());
189 
190  // gather operation for each cell in workset
191  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
192  std::vector<GlobalOrdinal> GIDs;
193  std::vector<int> LIDs;
194  std::size_t cellLocalId = localCellIds[worksetCellIndex];
195 
196  dofMngr.getElementGIDs(cellLocalId,GIDs);
197 
198  // caculate the local IDs for this element
199  LIDs.resize(GIDs.size());
200  for(std::size_t i=0;i<GIDs.size();i++)
201  LIDs[i] = x.Map().LID(GIDs[i]);
202 
203  // loop over basis functions and fill the fields
204  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
205  int offset = elmtOffset[basis];
206  int lid = LIDs[offset];
207  fc[fieldStr](worksetCellIndex,basis) = x[lid];
208  }
209  }
210  }
211 }
212 
213 #ifdef PANZER_HAVE_FEI
214 void scatter_to_vector(const std::string & blockId, const panzer::DOFManagerFEI<int,int> & dofMngr,
215  const std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc,
216  const std::vector<std::size_t> & localCellIds,
217  Epetra_Vector & x)
218 {
219 
220  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> >::const_iterator fieldItr;
221  for(fieldItr=fc.begin();fieldItr!=fc.end();++fieldItr) {
222  std::string fieldStr = fieldItr->first;
223  int fieldNum = dofMngr.getFieldNum(fieldStr);
224  const Kokkos::DynRankView<double,PHX::Device> & data = fieldItr->second;
225 
226  // gather operation for each cell in workset
227  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
228  std::vector<int> GIDs, LIDs;
229  std::size_t cellLocalId = localCellIds[worksetCellIndex];
230 
231  dofMngr.getElementGIDs(cellLocalId,GIDs);
232 
233  // caculate the local IDs for this element
234  LIDs.resize(GIDs.size());
235  for(std::size_t i=0;i<GIDs.size();i++)
236  LIDs[i] = x.Map().LID(GIDs[i]);
237 
238  const std::vector<int> & elmtOffset = dofMngr.getGIDFieldOffsets(blockId,fieldNum);
239 
240  // loop over basis functions and fill the fields
241  for(int basis=0;basis<data.dimension(1);basis++) {
242  int offset = elmtOffset[basis];
243  int lid = LIDs[offset];
244  x[lid] = data(worksetCellIndex,basis);
245  }
246  }
247  }
248 }
249 #endif
250 
252  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > & localIds)
253 {
254  // defines ordering of blocks
255  std::vector<std::string> blockIds;
256  mesh.getElementBlockNames(blockIds);
257 
258  std::vector<std::string>::const_iterator idItr;
259  for(idItr=blockIds.begin();idItr!=blockIds.end();++idItr) {
260  std::string blockId = *idItr;
261 
262  localIds[blockId] = Teuchos::rcp(new std::vector<std::size_t>);
263  std::vector<std::size_t> & localBlockIds = *localIds[blockId];
264 
265  // grab elements on this block
266  std::vector<stk::mesh::Entity> blockElmts;
267  mesh.getMyElements(blockId,blockElmts);
268 
269  std::vector<stk::mesh::Entity>::const_iterator itr;
270  for(itr=blockElmts.begin();itr!=blockElmts.end();++itr)
271  localBlockIds.push_back(mesh.elementLocalId(*itr));
272 
273  std::sort(localBlockIds.begin(),localBlockIds.end());
274  }
275 }
276 
277 }
virtual void getElementGIDs(LocalOrdinalT localElmtId, std::vector< GlobalOrdinalT > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
void setSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
std::size_t elementLocalId(stk::mesh::Entity elmt) const
stk::mesh::Field< double > SolutionFieldType
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
void write_solution_data(const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &dofMngr, panzer_stk::STK_Interface &mesh, const Epetra_MultiVector &x, const std::string &prefix, const std::string &postfix)
virtual const std::vector< int > & getGIDFieldOffsets(const std::string &blockId, int fieldNum) const =0
Use the field pattern so that you can find a particular field in the GIDs array.
void getElementBlockNames(std::vector< std::string > &names) const
static void gather_in_block(const std::string &blockId, const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &dofMngr, const Epetra_Vector &x, const std::vector< std::size_t > &localCellIds, std::map< std::string, Kokkos::DynRankView< double, PHX::Device > > &fc)
template void write_solution_data< int >(const panzer::UniqueGlobalIndexer< int, int > &dofMngr, panzer_stk::STK_Interface &mesh, const Epetra_MultiVector &x, const std::string &prefix, const std::string &postfix)
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
static void build_local_ids(const panzer_stk::STK_Interface &mesh, std::map< std::string, Teuchos::RCP< std::vector< std::size_t > > > &localIds)
PHX::MDField< const ScalarT, Cell, IP > field
virtual const std::string & getFieldString(int num) const =0
Reverse lookup of the field string from a field number.
void getSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, ArrayT &solutionValues) const
void write_cell_data(panzer_stk::STK_Interface &mesh, const std::vector< double > &data, const std::string &fieldName)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
virtual const std::vector< int > & getBlockFieldNumbers(const std::string &blockId) const =0