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 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #include "PanzerAdaptersSTK_config.hpp"
12 
13 #ifdef PANZER_HAVE_EPETRA_STACK
14 
15 #include "Panzer_STK_Utilities.hpp"
16 #include "Panzer_GlobalIndexer.hpp"
17 
18 #include "Kokkos_DynRankView.hpp"
19 
20 #include <stk_mesh/base/FieldBase.hpp>
21 
22 namespace panzer_stk {
23 
24 static void gather_in_block(const std::string & blockId, const panzer::GlobalIndexer& dofMngr,
25  const Epetra_Vector & x,const std::vector<std::size_t> & localCellIds,
26  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc);
27 
28 static void build_local_ids(const panzer_stk::STK_Interface & mesh,
29  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > & localIds);
30 
31 void write_cell_data(panzer_stk::STK_Interface & mesh,const std::vector<double> & data,const std::string & fieldName)
32 {
33  std::vector<std::string> blocks;
34  mesh.getElementBlockNames(blocks);
35 
36  // loop over element blocks
37  for(std::size_t eb=0;eb<blocks.size();eb++) {
38  const std::string & blockId = blocks[eb];
40 
41  std::vector<stk::mesh::Entity> elements;
42  mesh.getMyElements(blockId,elements);
43 
44  // loop over elements in this block
45  for(std::size_t el=0;el<elements.size();el++) {
46  std::size_t localId = mesh.elementLocalId(elements[el]);
47  double * solnData = stk::mesh::field_data(*field,elements[el]);
48  TEUCHOS_ASSERT(solnData!=0); // sanity check
49  solnData[0] = data[localId];
50  }
51  }
52 }
53 
54 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)
55 {
56  write_solution_data(dofMngr,mesh,*x(0),prefix,postfix);
57 }
58 
59 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)
60 {
61  typedef Kokkos::DynRankView<double,PHX::Device> FieldContainer;
62 
63  // get local IDs
64  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > localIds;
65  build_local_ids(mesh,localIds);
66 
67  // loop over all element blocks
68  for(const auto & itr : localIds) {
69  const auto blockId = itr.first;
70  const auto & localCellIds = *(itr.second);
71 
72  std::map<std::string,FieldContainer> data;
73 
74  // get all solution data for this block
75  gather_in_block(blockId,dofMngr,x,localCellIds,data);
76 
77  // write out to stk mesh
78  for(const auto & dataItr : data)
79  mesh.setSolutionFieldData(prefix+dataItr.first+postfix,blockId,localCellIds,dataItr.second);
80  }
81 }
82 
83 void gather_in_block(const std::string & blockId, const panzer::GlobalIndexer& dofMngr,
84  const Epetra_Vector & x,const std::vector<std::size_t> & localCellIds,
85  std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > & fc)
86 {
87  const std::vector<int> & fieldNums = dofMngr.getBlockFieldNumbers(blockId);
88 
89  for(std::size_t fieldIndex=0;fieldIndex<fieldNums.size();fieldIndex++) {
90  int fieldNum = fieldNums[fieldIndex];
91  std::string fieldStr = dofMngr.getFieldString(fieldNum);
92 
93  // grab the field
94  const std::vector<int> & elmtOffset = dofMngr.getGIDFieldOffsets(blockId,fieldNum);
95  fc[fieldStr] = Kokkos::DynRankView<double,PHX::Device>("fc",localCellIds.size(),elmtOffset.size());
96  auto field = Kokkos::create_mirror_view(fc[fieldStr]);
97 
98 
99  // gather operation for each cell in workset
100  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
101  std::vector<panzer::GlobalOrdinal> GIDs;
102  std::vector<int> LIDs;
103  std::size_t cellLocalId = localCellIds[worksetCellIndex];
104 
105  dofMngr.getElementGIDs(cellLocalId,GIDs);
106 
107  // caculate the local IDs for this element
108  LIDs.resize(GIDs.size());
109  for(std::size_t i=0;i<GIDs.size();i++)
110  LIDs[i] = x.Map().LID(GIDs[i]);
111 
112  // loop over basis functions and fill the fields
113  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
114  int offset = elmtOffset[basis];
115  int lid = LIDs[offset];
116  field(worksetCellIndex,basis) = x[lid];
117  }
118  }
119  Kokkos::deep_copy(fc[fieldStr], field);
120  }
121 }
122 
123 void build_local_ids(const panzer_stk::STK_Interface & mesh,
124  std::map<std::string,Teuchos::RCP<std::vector<std::size_t> > > & localIds)
125 {
126  // defines ordering of blocks
127  std::vector<std::string> blockIds;
128  mesh.getElementBlockNames(blockIds);
129 
130  std::vector<std::string>::const_iterator idItr;
131  for(idItr=blockIds.begin();idItr!=blockIds.end();++idItr) {
132  std::string blockId = *idItr;
133 
134  localIds[blockId] = Teuchos::rcp(new std::vector<std::size_t>);
135  std::vector<std::size_t> & localBlockIds = *localIds[blockId];
136 
137  // grab elements on this block
138  std::vector<stk::mesh::Entity> blockElmts;
139  mesh.getMyElements(blockId,blockElmts);
140 
141  std::vector<stk::mesh::Entity>::const_iterator itr;
142  for(itr=blockElmts.begin();itr!=blockElmts.end();++itr)
143  localBlockIds.push_back(mesh.elementLocalId(*itr));
144 
145  std::sort(localBlockIds.begin(),localBlockIds.end());
146  }
147 }
148 
149 }
150 
151 #endif // PANZER_HAVE_EPETRA_STACK
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)
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