Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_ParameterListCallbackBlocked.hpp
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 #ifndef __Panzer_STK_ParameterListCallbackBlocked_hpp__
12 #define __Panzer_STK_ParameterListCallbackBlocked_hpp__
13 
14 #include "PanzerAdaptersSTK_config.hpp"
15 #ifdef PANZER_HAVE_TEKO
16 
17 #include "Teuchos_RCP.hpp"
19 
20 #include "Teko_RequestCallback.hpp"
21 
24 #ifdef PANZER_HAVE_EPETRA_STACK
26 #endif
28 
29 #include <vector>
30 #include <map>
31 
32 namespace panzer_stk {
33 
34 class STKConnManager;
35 
40 class ParameterListCallbackBlocked : public Teko::RequestCallback<Teuchos::RCP<Teuchos::ParameterList> > {
41 public:
42  ParameterListCallbackBlocked(const Teuchos::RCP<const panzer_stk::STKConnManager> & connManager,
44  const Teuchos::RCP<const panzer::BlockedDOFManager> & auxBlkDofs=Teuchos::null);
45 
46  Teuchos::RCP<Teuchos::ParameterList> request(const Teko::RequestMesg & rm);
47 
48  bool handlesRequest(const Teko::RequestMesg & rm);
49 
50  void preRequest(const Teko::RequestMesg & rm);
51 
52 private:
53 
54  bool isField(const std::string& field) const
55  {
56  // Check both the main and auxiliary UGIs.
57  bool useAux(true);
58  std::vector<Teuchos::RCP<panzer::GlobalIndexer>>
59  fieldDOFMngrs = blocked_ugi_->getFieldDOFManagers();
60  for (int b(0); b < static_cast<int>(fieldDOFMngrs.size()); ++b)
61  {
62  for (int f(0); f < fieldDOFMngrs[b]->getNumFields(); ++f)
63  {
64  if (fieldDOFMngrs[b]->getFieldString(f) == field)
65  useAux = false;
66  }
67  }
68  if (useAux)
69  return (aux_blocked_ugi_->getFieldNum(field) != -1);
70  else
71  return (blocked_ugi_->getFieldNum(field) != -1);
72  }
73 
74  void buildArrayToVectorTpetra(int block,const std::string & field, const bool useAux = false);
75  void buildCoordinatesTpetra(const std::string & field, const bool useAux = false);
76 
77 #ifdef PANZER_HAVE_EPETRA_STACK
78  void buildArrayToVectorEpetra(int block,const std::string & field, const bool useAux = false);
79  void buildCoordinatesEpetra(const std::string & field, const bool useAux = false);
80 #endif
81 
82  // this method assumes handlesRequest(rm)==true
83  std::string getHandledField(const Teuchos::ParameterList & pl) const;
84 
85  void setFieldByKey(const std::string & key,const std::string & field,Teuchos::ParameterList & pl) const;
86 
87  // Get the coordinate vector by field, note that it will check to make sure "buildCoordinates" has
88  // been called.
89  const std::vector<double> & getCoordinateByField(int dim,const std::string & field) const;
90 
91  // Access the field pattern associated with this field (this will not work if the field pattern
92  // and element blocks are different
93  Teuchos::RCP<const panzer::Intrepid2FieldPattern> getFieldPattern(const std::string & fieldName, const bool useAux = false) const;
94 
95  // Generally used members
99 
100  std::map<std::string,Teuchos::RCP<const panzer::Intrepid2FieldPattern> > fieldPatterns_;
101 
102  // look up by field name (field name to coordinates
103 
104  std::map<std::string,std::vector<double> > xcoords_;
105  std::map<std::string,std::vector<double> > ycoords_;
106  std::map<std::string,std::vector<double> > zcoords_;
107 
108  mutable std::map<std::string,Teuchos::RCP<const panzer::ArrayToFieldVector> > arrayToVectorTpetra_;
109 #ifdef PANZER_HAVE_EPETRA_STACK
110  mutable std::map<std::string,Teuchos::RCP<const panzer::ArrayToFieldVectorEpetra> > arrayToVectorEpetra_;
111 #endif
112 
114 #ifdef PANZER_HAVE_EPETRA_STACK
116 #endif
117 
118  bool returnTpetraObjects_;
119 };
120 
121 }
122 
123 #endif // PANZER_HAVE_TEKO
124 
125 #endif
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...