Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_TpetraLinearObjFactory_decl.hpp
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 #ifndef __Panzer_TpetraLinearObjFactory_decl_hpp__
44 #define __Panzer_TpetraLinearObjFactory_decl_hpp__
45 
46 #include <map>
47 
48 // Tpetra includes
49 #include "Tpetra_Map.hpp"
50 #include "Tpetra_CrsGraph.hpp"
51 #include "Tpetra_CrsMatrix.hpp"
52 #include "Tpetra_Import.hpp"
53 #include "Tpetra_Export.hpp"
54 
55 #include "PanzerDiscFE_config.hpp"
56 #include "Panzer_GlobalIndexer.hpp"
59 #include "Panzer_ScatterResidual_Tpetra.hpp"
60 #include "Panzer_ScatterDirichletResidual_Tpetra.hpp"
61 #include "Panzer_GatherSolution_Tpetra.hpp"
62 #include "Panzer_GatherTangent_Tpetra.hpp"
63 #include "Panzer_GatherOrientation.hpp"
66 
67 #include"Panzer_NodeType.hpp"
68 
69 #include "Teuchos_RCP.hpp"
71 
72 namespace panzer {
73 
74 template <typename Traits,typename ScalarT,typename LocalOrdinalT,typename GlobalOrdinalT,typename NodeT=panzer::TpetraNodeType>
76  , public ThyraObjFactory<ScalarT> {
77 public:
79  typedef Tpetra::Vector<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> VectorType;
80  typedef Tpetra::CrsMatrix<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> CrsMatrixType;
81  typedef Tpetra::CrsGraph<LocalOrdinalT,GlobalOrdinalT,NodeT> CrsGraphType;
82  typedef Tpetra::Map<LocalOrdinalT,GlobalOrdinalT,NodeT> MapType;
83  typedef Tpetra::Import<LocalOrdinalT,GlobalOrdinalT,NodeT> ImportType;
84  typedef Tpetra::Export<LocalOrdinalT,GlobalOrdinalT,NodeT> ExportType;
85 
87  const Teuchos::RCP<const GlobalIndexer> & gidProvider);
88 
90  const Teuchos::RCP<const GlobalIndexer> & rowProvider,
91  const Teuchos::RCP<const GlobalIndexer> & colProvider);
92 
93  virtual ~TpetraLinearObjFactory();
94 
95 /*************** Linear object factory methods *******************/
96 
97  virtual void readVector(const std::string & /* identifier */, LinearObjContainer & /* loc */, int /* id */) const
98  { TEUCHOS_ASSERT(false); }
99 
100  virtual void writeVector(const std::string & /* identifier */, const LinearObjContainer & /* loc */, int /* id */) const
101  { TEUCHOS_ASSERT(false); }
102 
104 
106  { return buildLinearObjContainer(); }
107 
109 
111  { return buildGhostedLinearObjContainer(); }
112 
113  virtual void globalToGhostContainer(const LinearObjContainer & container,
114  LinearObjContainer & ghostContainer,int) const;
115  virtual void ghostToGlobalContainer(const LinearObjContainer & ghostContainer,
116  LinearObjContainer & container,int) const;
117 
124  virtual void adjustForDirichletConditions(const LinearObjContainer & localBCRows,
125  const LinearObjContainer & globalBCRows,
126  LinearObjContainer & ghostedObjs,
127  bool zeroVectorRows=false, bool adjustX=false) const;
128 
132  virtual void applyDirichletBCs(const LinearObjContainer & counter,
133  LinearObjContainer & result) const;
134 
140 
141 #ifdef PANZER_HAVE_EPETRA_STACK
142 
146  virtual Teuchos::RCP<WriteVector_GlobalEvaluationData> buildWriteDomainContainer() const;
147 #endif
148 
151  virtual Teuchos::MpiComm<int> getComm() const;
152 
154  template <typename EvalT>
157 
159  template <typename EvalT>
162 
164  template <typename EvalT>
167 
169  template <typename EvalT>
172 
174  template <typename EvalT>
177 
179  template <typename EvalT>
182 
183 /*************** From ThyraObjFactory *******************/
184 
187 
190 
193 
194 /*************** Tpetra construction functions *******************/
195 
202 
203 /*************** Generic helper functions for container setup *******************/
204 
210  void initializeContainer(int,LinearObjContainer & loc) const;
211 
218 
224  void initializeGhostedContainer(int,LinearObjContainer & loc) const;
225 
232 
233 /*************** Tpetra based methods *******************/
234 
238 
242 
245 
248 
252 
256 
259 
262  { return gidProvider_; }
263 
266  { return gidProvider_; }
267 
268  virtual void beginFill(LinearObjContainer & loc) const;
269  virtual void endFill(LinearObjContainer & loc) const;
270 
271 protected:
272 
273  void ghostToGlobalTpetraVector(const Tpetra::Vector<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> & in,
274  Tpetra::Vector<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> & out, bool col) const;
275  void ghostToGlobalTpetraMatrix(const Tpetra::CrsMatrix<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> & in,
276  Tpetra::CrsMatrix<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> & out) const;
277  void globalToGhostTpetraVector(const Tpetra::Vector<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT>& in,
278  Tpetra::Vector<ScalarT,LocalOrdinalT,GlobalOrdinalT,NodeT> & out, bool col) const;
279 
280  // get the map from the matrix
285 
286  // get the graph of the crs matrix
289 
290  // storage for Tpetra graphs and maps
302 
305 
307 
310 };
311 
312 }
313 
314 #endif
Teuchos::RCP< panzer::CloneableEvaluator > buildGatherOrientation() const
Use preconstructed gather evaluators.
virtual const Teuchos::RCP< Tpetra::Export< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedExport() const
get exporter for converting an overalapped object to a &quot;normal&quot; object
virtual void globalToGhostContainer(const LinearObjContainer &container, LinearObjContainer &ghostContainer, int) const
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > domainSpace_
Teuchos::RCP< Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > > getTpetraVector() const
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedColMap() const
Teuchos::RCP< const GlobalIndexer > gidProvider_
virtual Teuchos::MpiComm< int > getComm() const
Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > cGhostedMap_
virtual Teuchos::RCP< LinearObjContainer > buildGhostedLinearObjContainer() const
void initializeGhostedContainer(int, LinearObjContainer &loc) const
virtual Teuchos::RCP< ReadOnlyVector_GlobalEvaluationData > buildReadOnlyDomainContainer() const
Teuchos::RCP< panzer::CloneableEvaluator > buildGatherDomain() const
Use preconstructed gather evaluators.
Teuchos::RCP< Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedTpetraVector() const
virtual const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > > buildGhostedGraph() const
virtual Teuchos::RCP< const Thyra::VectorSpaceBase< ScalarT > > getThyraDomainSpace() const
Get the domain space.
Teuchos::RCP< Tpetra::CrsMatrix< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > > getTpetraMatrix() const
Teuchos::RCP< panzer::CloneableEvaluator > buildGatherTangent() const
Use preconstructed gather evaluators.
virtual const Teuchos::RCP< Tpetra::Import< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedImport() const
get importer for converting an overalapped object to a &quot;normal&quot; object
Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > cMap_
void initializeContainer(int, LinearObjContainer &loc) const
Teuchos::RCP< Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > > getTpetraColVector() const
Teuchos::RCP< const panzer::GlobalIndexer > getDomainGlobalIndexer() const
Get the domain global indexer this factory was created with.
virtual void endFill(LinearObjContainer &loc) const
Gathers tangent vectors dx/dp for computing df/dx*dx/dp + df/dp into the nodal fields of the field ma...
virtual const Teuchos::RCP< Tpetra::Import< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedColImport() const
virtual void beginFill(LinearObjContainer &loc) const
Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > MapType
virtual const Teuchos::RCP< const Teuchos::Comm< int > > getTeuchosComm() const
get exporter for converting an overalapped object to a &quot;normal&quot; object
Tpetra::Import< LocalOrdinalT, GlobalOrdinalT, NodeT > ImportType
virtual const Teuchos::RCP< Tpetra::Export< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedColExport() const
Teuchos::RCP< Tpetra::CrsMatrix< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedTpetraMatrix() const
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > getColMap() const
virtual void adjustForDirichletConditions(const LinearObjContainer &localBCRows, const LinearObjContainer &globalBCRows, LinearObjContainer &ghostedObjs, bool zeroVectorRows=false, bool adjustX=false) const
void globalToGhostTpetraVector(const Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > &in, Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > &out, bool col) const
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we&#39;re performing.
void ghostToGlobalTpetraMatrix(const Tpetra::CrsMatrix< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > &in, Tpetra::CrsMatrix< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > &out) const
Teuchos::RCP< panzer::CloneableEvaluator > buildGather() const
Use preconstructed gather evaluators.
Tpetra::Export< LocalOrdinalT, GlobalOrdinalT, NodeT > ExportType
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Tpetra::CrsMatrix< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > CrsMatrixType
virtual Teuchos::RCP< LinearObjContainer > buildLinearObjContainer() const
Teuchos::RCP< panzer::CloneableEvaluator > buildScatterDirichlet() const
Use preconstructed dirichlet scatter evaluators.
TpetraLinearObjContainer< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > ContainerType
Gathers solution values from the Newton solution vector into the nodal fields of the field manager...
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedMap() const
get the ghosted map from the matrix
Pushes residual values into the residual vector for a Newton-based solve.
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > buildColMap() const
virtual const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > > buildGraph() const
Teuchos::RCP< Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedTpetraColVector() const
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > getMap() const
get the map from the matrix
Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > ghostedMap_
Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > > graph_
virtual const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGraph() const
get the graph of the crs matrix
Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > > ghostedGraph_
Pushes residual values into the residual vector for a Newton-based solve.
TpetraLinearObjFactory(const Teuchos::RCP< const Teuchos::Comm< int > > &comm, const Teuchos::RCP< const GlobalIndexer > &gidProvider)
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > buildMap() const
Teuchos::RCP< const GlobalIndexer > colGidProvider_
virtual Teuchos::RCP< LinearObjContainer > buildPrimitiveGhostedLinearObjContainer() const
Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > map_
virtual Teuchos::RCP< const Thyra::VectorSpaceBase< ScalarT > > getThyraRangeSpace() const
Get the range space.
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > buildGhostedColMap() const
virtual void ghostToGlobalContainer(const LinearObjContainer &ghostContainer, LinearObjContainer &container, int) const
Gathers orientations per field from the global indexer and stores them in the field manager...
virtual void readVector(const std::string &, LinearObjContainer &, int) const
virtual Teuchos::RCP< Thyra::LinearOpBase< ScalarT > > getThyraMatrix() const
Get a matrix operator.
Teuchos::RCP< const Teuchos::Comm< int > > comm_
#define TEUCHOS_ASSERT(assertion_test)
void ghostToGlobalTpetraVector(const Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > &in, Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > &out, bool col) const
Teuchos::RCP< panzer::CloneableEvaluator > buildScatter() const
Use preconstructed scatter evaluators.
virtual void writeVector(const std::string &, const LinearObjContainer &, int) const
virtual Teuchos::RCP< LinearObjContainer > buildPrimitiveLinearObjContainer() const
Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > CrsGraphType
Tpetra::Vector< ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT > VectorType
Teuchos::RCP< const panzer::GlobalIndexer > getRangeGlobalIndexer() const
Get the domain global indexer this factory was created with.
virtual const Teuchos::RCP< Tpetra::CrsGraph< LocalOrdinalT, GlobalOrdinalT, NodeT > > getGhostedGraph() const
get the ghosted graph of the crs matrix
virtual const Teuchos::RCP< Tpetra::Map< LocalOrdinalT, GlobalOrdinalT, NodeT > > buildGhostedMap() const
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > rangeSpace_
virtual void applyDirichletBCs(const LinearObjContainer &counter, LinearObjContainer &result) const