MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_BlockedCoordinatesTransferFactory_def.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_BLOCKEDCOORDINATESTRANSFER_FACTORY_DEF_HPP
47#define MUELU_BLOCKEDCOORDINATESTRANSFER_FACTORY_DEF_HPP
48
49#include "Xpetra_ImportFactory.hpp"
50#include "Xpetra_MultiVectorFactory.hpp"
51#include "Xpetra_MapFactory.hpp"
52#include "Xpetra_IO.hpp"
53
55
56#include "MueLu_Level.hpp"
57#include "MueLu_Monitor.hpp"
58
59namespace MueLu {
60
61 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63 RCP<ParameterList> validParamList = rcp(new ParameterList());
64
65 validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Factory for coordinates generation");
66 validParamList->set<RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
67 return validParamList;
68 }
69
70 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72 Input(coarseLevel, "CoarseMap");
73
74 // Make sure the Level knows I need these sub-Factories
75 const size_t numSubFactories = NumFactories();
76 for(size_t i=0; i<numSubFactories; i++) {
77 const RCP<const FactoryBase>& myFactory = subFactories_[i];
78 coarseLevel.DeclareInput("Coordinates", myFactory.getRawPtr(), this);
79 }
80
81 // call DeclareInput of all user-given transfer factories
82 for (std::vector<RCP<const FactoryBase> >::const_iterator it = subFactories_.begin(); it != subFactories_.end(); ++it)
83 (*it)->CallDeclareInput(coarseLevel);
84 }
85
86 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
88 FactoryMonitor m(*this, "Build", coarseLevel);
89
90 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO> dMV;
91 typedef Xpetra::BlockedMultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO> dBV;
92
93 GetOStream(Runtime0) << "Transferring (blocked) coordinates" << std::endl;
94
95 const size_t numSubFactories = NumFactories();
96 std::vector<RCP<const Map> > subBlockMaps(numSubFactories);
97 std::vector<RCP<dMV> > subBlockCoords(numSubFactories);
98
99 if (coarseLevel.IsAvailable("Coordinates", this)) {
100 GetOStream(Runtime0) << "Reusing coordinates" << std::endl;
101 return;
102 }
103
104 // Get components
105 for(size_t i=0; i<numSubFactories; i++) {
106 GetOStream(Runtime1) << "Generating Coordinates for block " << i <<"/"<<numSubFactories <<std::endl;
107 const RCP<const FactoryBase>& myFactory = subFactories_[i];
108 myFactory->CallBuild(coarseLevel);
109 subBlockCoords[i] = coarseLevel.Get<RCP<dMV> >("Coordinates", myFactory.get());
110 subBlockMaps[i] = subBlockCoords[i]->getMap();
111 }
112
113 // Blocked Map
114 RCP<const BlockedMap> coarseCoordMapBlocked;
115
116 {
117 // coarseMap is being used to set up the domain map of tentative P, and therefore, the row map of Ac
118 // Therefore, if we amalgamate coarseMap, logical nodes in the coordinates vector would correspond to
119 // logical blocks in the matrix
120 RCP<const BlockedMap> coarseMap = Get< RCP<const BlockedMap> >(coarseLevel, "CoarseMap");
121 bool thyraMode = coarseMap->getThyraMode();
122
123 ArrayView<const GO> elementAList = coarseMap->getFullMap()->getLocalElementList();
124
125 LO blkSize = 1;
126 if (rcp_dynamic_cast<const StridedMap>(coarseMap->getMap(0, thyraMode)) != Teuchos::null)
127 blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap->getMap(0, thyraMode))->getFixedBlockSize();
128
129 for(size_t i=1; i<numSubFactories; i++) {
130 LO otherBlkSize = 1;
131 if (rcp_dynamic_cast<const StridedMap>(coarseMap->getMap(i, thyraMode)) != Teuchos::null)
132 otherBlkSize = rcp_dynamic_cast<const StridedMap>(coarseMap->getMap(i, thyraMode))->getFixedBlockSize();
133 TEUCHOS_TEST_FOR_EXCEPTION(otherBlkSize != blkSize, Exceptions::RuntimeError, "BlockedCoordinatesTransferFactory: Subblocks have different Block sizes. This is not yet supported.");
134 }
135
136 GO indexBase = coarseMap->getFullMap()->getIndexBase();
137 size_t numElements = elementAList.size() / blkSize;
138 Array<GO> elementList(numElements);
139
140 // Amalgamate the map
141 for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
142 elementList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
143
144 RCP<const Map> coarseCoordMap = MapFactory::Build(coarseMap->getFullMap()->lib(),
145 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), elementList, indexBase, coarseMap->getFullMap()->getComm());
146
147 coarseCoordMapBlocked = rcp(new BlockedMap(coarseCoordMap, subBlockMaps, thyraMode));
148 }
149
150 // Build blocked coordinates vector
151 RCP<dBV> bcoarseCoords = rcp(new dBV(coarseCoordMapBlocked,subBlockCoords));
152
153 // Turn the blocked coordinates vector into an unblocked one
154 RCP<dMV> coarseCoords = bcoarseCoords->Merge();
155 Set<RCP<dMV> >(coarseLevel, "Coordinates", coarseCoords);
156 }
157
158 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
160 subFactories_.push_back(factory);
161 }
162
163
164
165} // namespace MueLu
166
167#endif // MUELU_BLOCKEDCOORDINATESTRANSFER_FACTORY_DEF_HPP
std::vector< RCP< const FactoryBase > > subFactories_
list of user-defined sub Factories
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
void AddFactory(const RCP< const FactoryBase > &factory)
Add (sub) coords factory in the end of list of factories in BlockedCoordinatesTransferFactory.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput().
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose).