Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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"
46 #include "Panzer_GlobalIndexer.hpp"
47 
48 #include "Kokkos_DynRankView.hpp"
49 
50 #include <stk_mesh/base/FieldBase.hpp>
51 
52 namespace panzer_stk {
53 
54 static void gather_in_block(const std::string & blockId, const panzer::GlobalIndexer& dofMngr,
55  const Epetra_Vector & x,const std::vector<std::size_t> & localCellIds,
56  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc);
57 
58 static void build_local_ids(const panzer_stk::STK_Interface & mesh,
59  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > & localIds);
60 
61 void write_cell_data(panzer_stk::STK_Interface & mesh,const std::vector<double> & data,const std::string & fieldName)
62 {
63  std::vector<std::string> blocks;
64  mesh.getElementBlockNames(blocks);
65 
66  // loop over element blocks
67  for(std::size_t eb=0;eb<blocks.size();eb++) {
68  const std::string & blockId = blocks[eb];
70 
71  std::vector<stk::mesh::Entity> elements;
72  mesh.getMyElements(blockId,elements);
73 
74  // loop over elements in this block
75  for(std::size_t el=0;el<elements.size();el++) {
76  std::size_t localId = mesh.elementLocalId(elements[el]);
77  double * solnData = stk::mesh::field_data(*field,elements[el]);
78  TEUCHOS_ASSERT(solnData!=0); // sanity check
79  solnData[0] = data[localId];
80  }
81  }
82 }
83 
84 void write_solution_data(const panzer::GlobalIndexer& dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_MultiVector & x,const std::string & prefix,const std::string & postfix)
85 {
86  write_solution_data(dofMngr,mesh,*x(0),prefix,postfix);
87 }
88 
89 void write_solution_data(const panzer::GlobalIndexer& dofMngr,panzer_stk::STK_Interface & mesh,const Epetra_Vector & x,const std::string & prefix,const std::string & postfix)
90 {
91  typedef Kokkos::DynRankView<double,PHX::Device> FieldContainer;
92 
93  // get local IDs
94  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > localIds;
95  build_local_ids(mesh,localIds);
96 
97  // loop over all element blocks
98  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > >::const_iterator itr;
99  for(itr=localIds.begin();itr!=localIds.end();++itr) {
100  std::string blockId = itr->first;
101  const std::vector<std::size_t> & localCellIds = *(itr->second);
102 
103  std::map<std::string,FieldContainer> data;
104 
105  // get all solution data for this block
106  gather_in_block(blockId,dofMngr,x,localCellIds,data);
107 
108  // write out to stk mesh
109  std::map<std::string,FieldContainer>::iterator dataItr;
110  for(dataItr=data.begin();dataItr!=data.end();++dataItr)
111  mesh.setSolutionFieldData(prefix+dataItr->first+postfix,blockId,localCellIds,dataItr->second);
112  }
113 }
114 
115 void gather_in_block(const std::string & blockId, const panzer::GlobalIndexer& dofMngr,
116  const Epetra_Vector & x,const std::vector<std::size_t> & localCellIds,
117  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc)
118 {
119  const std::vector<int> & fieldNums = dofMngr.getBlockFieldNumbers(blockId);
120 
121  for(std::size_t fieldIndex=0;fieldIndex<fieldNums.size();fieldIndex++) {
122  int fieldNum = fieldNums[fieldIndex];
123  std::string fieldStr = dofMngr.getFieldString(fieldNum);
124 
125  // grab the field
126  const std::vector<int> & elmtOffset = dofMngr.getGIDFieldOffsets(blockId,fieldNum);
127  fc[fieldStr] = Kokkos::DynRankView<double,PHX::Device>("fc",localCellIds.size(),elmtOffset.size());
128 
129  // gather operation for each cell in workset
130  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
131  std::vector<panzer::GlobalOrdinal> GIDs;
132  std::vector<int> LIDs;
133  std::size_t cellLocalId = localCellIds[worksetCellIndex];
134 
135  dofMngr.getElementGIDs(cellLocalId,GIDs);
136 
137  // caculate the local IDs for this element
138  LIDs.resize(GIDs.size());
139  for(std::size_t i=0;i<GIDs.size();i++)
140  LIDs[i] = x.Map().LID(GIDs[i]);
141 
142  // loop over basis functions and fill the fields
143  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
144  int offset = elmtOffset[basis];
145  int lid = LIDs[offset];
146  fc[fieldStr](worksetCellIndex,basis) = x[lid];
147  }
148  }
149  }
150 }
151 
153  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > & localIds)
154 {
155  // defines ordering of blocks
156  std::vector<std::string> blockIds;
157  mesh.getElementBlockNames(blockIds);
158 
159  std::vector<std::string>::const_iterator idItr;
160  for(idItr=blockIds.begin();idItr!=blockIds.end();++idItr) {
161  std::string blockId = *idItr;
162 
163  localIds[blockId] = Teuchos::rcp(new std::vector<std::size_t>);
164  std::vector<std::size_t> & localBlockIds = *localIds[blockId];
165 
166  // grab elements on this block
167  std::vector<stk::mesh::Entity> blockElmts;
168  mesh.getMyElements(blockId,blockElmts);
169 
170  std::vector<stk::mesh::Entity>::const_iterator itr;
171  for(itr=blockElmts.begin();itr!=blockElmts.end();++itr)
172  localBlockIds.push_back(mesh.elementLocalId(*itr));
173 
174  std::sort(localBlockIds.begin(),localBlockIds.end());
175  }
176 }
177 
178 }
virtual const std::string & getFieldString(int num) const =0
Reverse lookup of the field string from a field number.
void setSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
void getElementBlockNames(std::vector< std::string > &names) const
stk::mesh::Field< double > SolutionFieldType
std::size_t elementLocalId(stk::mesh::Entity elmt) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void gather_in_block(const std::string &blockId, const panzer::GlobalIndexer &dofMngr, const Epetra_Vector &x, const std::vector< std::size_t > &localCellIds, std::map< std::string, Kokkos::DynRankView< double, PHX::Device > > &fc)
void write_solution_data(const panzer::GlobalIndexer &dofMngr, panzer_stk::STK_Interface &mesh, const Epetra_MultiVector &x, const std::string &prefix, const std::string &postfix)
static void build_local_ids(const panzer_stk::STK_Interface &mesh, std::map< std::string, Teuchos::RCP< std::vector< std::size_t > > > &localIds)
void write_cell_data(panzer_stk::STK_Interface &mesh, const std::vector< double > &data, const std::string &fieldName)
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
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.
virtual const std::vector< int > & getBlockFieldNumbers(const std::string &blockId) const =0
virtual void getElementGIDs(panzer::LocalOrdinal localElmtId, std::vector< panzer::GlobalOrdinal > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
#define TEUCHOS_ASSERT(assertion_test)
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const