MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_GeometricInterpolationPFactory_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 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
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_GEOMETRICINTERPOLATIONPFACTORY_DECL_HPP
47 #define MUELU_GEOMETRICINTERPOLATIONPFACTORY_DECL_HPP
48 
49 // Teuchos includes for dense linear algebra
53 
54 #include "Xpetra_CrsGraph.hpp"
55 
56 #include "MueLu_PFactory.hpp"
57 #include "MueLu_Level_fwd.hpp"
58 #include "MueLu_Aggregates_fwd.hpp"
59 
60 namespace MueLu{
61 
62  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64 #undef MUELU_GEOMETRICINTERPOLATIONPFACTORY_SHORT
65 #include "MueLu_UseShortNames.hpp"
66 
67  public:
68 
69  // Declare useful types
71  using realvaluedmultivector_type = Xpetra::MultiVector<real_type,LO,GO,Node>;
72 
74 
75 
78 
82 
84 
86 
87 
88  void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
89 
91 
93 
94 
95  void Build (Level& fineLevel, Level& coarseLevel) const;
96  void BuildP(Level& fineLevel, Level& coarseLevel) const;
97 
99 
100  private:
101  void BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A) const;
102  void BuildLinearP(Level& coarseLevel,
103  RCP<Matrix>& A, RCP<const CrsGraph>& prolongatorGraph,
104  RCP<realvaluedmultivector_type>& fineCoordinates,
105  RCP<realvaluedmultivector_type>& ghostCoordinates,
106  const int numDimensions, const bool keepD2,
107  RCP<Matrix>& P) const;
108  void ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints,
109  const Array<Array<real_type> > coord,
110  Array<real_type>& stencil) const;
111  void GetInterpolationFunctions(const LO numDimensions,
112  const Teuchos::SerialDenseVector<LO,real_type> parametricCoordinates,
113  real_type functions[4][8]) const;
114 
115  }; // class GeometricInterpolationPFactory
116 
117 } // namespace MueLu
118 
119 #define MUELU_GEOMETRICINTERPOLATIONPFACTORY_SHORT
120 #endif // MUELU_GEOMETRICINTERPOLATIONPFACTORY_DECL_HPP
void BuildLinearP(Level &coarseLevel, RCP< Matrix > &A, RCP< const CrsGraph > &prolongatorGraph, RCP< realvaluedmultivector_type > &fineCoordinates, RCP< realvaluedmultivector_type > &ghostCoordinates, const int numDimensions, const bool keepD2, RCP< Matrix > &P) const
void GetInterpolationFunctions(const LO numDimensions, const Teuchos::SerialDenseVector< LO, real_type > parametricCoordinates, real_type functions[4][8]) const
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints, const Array< Array< real_type > > coord, Array< real_type > &stencil) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Xpetra::MultiVector< real_type, LO, GO, Node > realvaluedmultivector_type
Factory that provides an interface for a concrete implementation of a prolongation operator...
void BuildConstantP(RCP< Matrix > &P, RCP< const CrsGraph > &prolongatorGraph, RCP< Matrix > &A) const
typename Teuchos::ScalarTraits< SC >::coordinateType real_type
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.