9 #ifndef MUELU_FACTORYFACTORY_DEF_HPP
10 #define MUELU_FACTORYFACTORY_DEF_HPP
13 #include "MueLu_AggregateQualityEstimateFactory.hpp"
14 #include "MueLu_AggregationExportFactory.hpp"
15 #include "MueLu_AmalgamationFactory.hpp"
16 #include "MueLu_BlackBoxPFactory.hpp"
17 #include "MueLu_BlockedCoarseMapFactory.hpp"
18 #include "MueLu_BlockedCoordinatesTransferFactory.hpp"
19 #include "MueLu_BlockedDirectSolver.hpp"
20 #include "MueLu_BlockedGaussSeidelSmoother.hpp"
21 #include "MueLu_BlockedJacobiSmoother.hpp"
22 #include "MueLu_BlockedPFactory.hpp"
23 #include "MueLu_BlockedRAPFactory.hpp"
24 #include "MueLu_BraessSarazinSmoother.hpp"
25 #include "MueLu_BrickAggregationFactory.hpp"
26 #include "MueLu_ClassicalMapFactory.hpp"
27 #include "MueLu_ClassicalPFactory.hpp"
28 #include "MueLu_CloneRepartitionInterface.hpp"
29 #include "MueLu_CoalesceDropFactory.hpp"
30 #include "MueLu_SmooVecCoalesceDropFactory.hpp"
31 #include "MueLu_CoarseMapFactory.hpp"
32 #include "MueLu_CoarseningVisualizationFactory.hpp"
33 #include "MueLu_ConstraintFactory.hpp"
34 #include "MueLu_CoordinatesTransferFactory.hpp"
35 #include "MueLu_DirectSolver.hpp"
36 #include "MueLu_DropNegativeEntriesFactory.hpp"
37 #include "MueLu_EminPFactory.hpp"
38 #include "MueLu_FilteredAFactory.hpp"
39 #include "MueLu_FineLevelInputDataFactory.hpp"
40 #include "MueLu_GeneralGeometricPFactory.hpp"
41 #include "MueLu_ReplicatePFactory.hpp"
42 #include "MueLu_CombinePFactory.hpp"
43 #include "MueLu_GenericRFactory.hpp"
44 #include "MueLu_GeometricInterpolationPFactory.hpp"
45 #include "MueLu_InterfaceAggregationFactory.hpp"
46 #include "MueLu_InterfaceMappingTransferFactory.hpp"
47 #include "MueLu_InitialBlockNumberFactory.hpp"
48 #include "MueLu_IndefBlockedDiagonalSmoother.hpp"
49 #include "MueLu_InverseApproximationFactory.hpp"
50 #include "MueLu_IsorropiaInterface.hpp"
51 #include "MueLu_LineDetectionFactory.hpp"
52 #include "MueLu_LocalOrdinalTransferFactory.hpp"
53 #include "MueLu_RepartitionInterface.hpp"
54 #include "MueLu_RepartitionBlockDiagonalFactory.hpp"
55 #include "MueLu_MapTransferFactory.hpp"
56 #include "MueLu_MatrixAnalysisFactory.hpp"
57 #include "MueLu_MultiVectorTransferFactory.hpp"
58 #include "MueLu_NotayAggregationFactory.hpp"
59 #include "MueLu_NullspaceFactory.hpp"
60 #include "MueLu_NullspacePresmoothFactory.hpp"
61 #include "MueLu_PatternFactory.hpp"
62 #include "MueLu_PgPFactory.hpp"
63 #include "MueLu_RebalanceBlockInterpolationFactory.hpp"
64 #include "MueLu_RebalanceBlockRestrictionFactory.hpp"
65 #include "MueLu_RebalanceBlockAcFactory.hpp"
66 #include "MueLu_RebalanceTransferFactory.hpp"
67 #include "MueLu_RegionRFactory.hpp"
68 #include "MueLu_RepartitionFactory.hpp"
69 #include "MueLu_RepartitionHeuristicFactory.hpp"
70 #include "MueLu_RAPFactory.hpp"
71 #include "MueLu_RAPShiftFactory.hpp"
72 #include "MueLu_RebalanceAcFactory.hpp"
73 #include "MueLu_ReorderBlockAFactory.hpp"
74 #include "MueLu_SaPFactory.hpp"
75 #include "MueLu_ScaledNullspaceFactory.hpp"
76 #include "MueLu_SegregatedAFactory.hpp"
77 #include "MueLu_SemiCoarsenPFactory.hpp"
78 #include "MueLu_SchurComplementFactory.hpp"
79 #include "MueLu_SimpleSmoother.hpp"
80 #include "MueLu_SmootherFactory.hpp"
81 #include "MueLu_StructuredAggregationFactory.hpp"
82 #include "MueLu_StructuredLineDetectionFactory.hpp"
83 #include "MueLu_SubBlockAFactory.hpp"
84 #ifdef HAVE_MUELU_TEKO
85 #include "MueLu_TekoSmoother.hpp"
87 #include "MueLu_TentativePFactory.hpp"
88 #include "MueLu_ToggleCoordinatesTransferFactory.hpp"
89 #include "MueLu_TogglePFactory.hpp"
90 #include "MueLu_TrilinosSmoother.hpp"
91 #include "MueLu_TransPFactory.hpp"
92 #include "MueLu_RfromP_Or_TransP.hpp"
93 #include "MueLu_UncoupledAggregationFactory.hpp"
94 #include "MueLu_HybridAggregationFactory.hpp"
95 #include "MueLu_UnsmooshFactory.hpp"
96 #include "MueLu_UserAggregationFactory.hpp"
97 #include "MueLu_UserPFactory.hpp"
98 #include "MueLu_UzawaSmoother.hpp"
99 #include "MueLu_VariableDofLaplacianFactory.hpp"
100 #include "MueLu_ZeroSubBlockAFactory.hpp"
101 #include "MueLu_ZoltanInterface.hpp"
102 #include "MueLu_Zoltan2Interface.hpp"
103 #include "MueLu_NodePartitionInterface.hpp"
105 #include "MueLu_CoalesceDropFactory_kokkos.hpp"
106 #include "MueLu_GeometricInterpolationPFactory_kokkos.hpp"
107 #ifdef HAVE_MUELU_DEPRECATED_CODE
108 #include "MueLu_SaPFactory_kokkos.hpp"
110 #include "MueLu_SemiCoarsenPFactory_kokkos.hpp"
111 #include "MueLu_StructuredAggregationFactory_kokkos.hpp"
112 #include "MueLu_TentativePFactory_kokkos.hpp"
113 #include "MueLu_MatrixFreeTentativePFactory.hpp"
114 #include "MueLu_RegionRFactory_kokkos.hpp"
116 #ifdef HAVE_MUELU_MATLAB
117 #include "MueLu_SingleLevelMatlabFactory.hpp"
118 #include "MueLu_TwoLevelMatlabFactory.hpp"
119 #include "MueLu_MatlabSmoother.hpp"
122 #ifdef HAVE_MUELU_INTREPID2
123 #include "MueLu_IntrepidPCoarsenFactory.hpp"
128 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
131 std::string factoryName;
134 factoryName = Teuchos::getValue<std::string>(param);
136 paramList = Teuchos::getValue<Teuchos::ParameterList>(param);
137 factoryName = paramList.
get<std::string>(
"factory");
141 if (factoryName ==
"AggregateQualityEstimateFactory")
return Build2<AggregateQualityEstimateFactory>(paramList, factoryMapIn, factoryManagersIn);
142 if (factoryName ==
"AggregationExportFactory")
return Build2<AggregationExportFactory>(paramList, factoryMapIn, factoryManagersIn);
143 if (factoryName ==
"AmalgamationFactory")
return Build2<AmalgamationFactory>(paramList, factoryMapIn, factoryManagersIn);
144 if (factoryName ==
"BlockedCoarseMapFactory")
return Build2<BlockedCoarseMapFactory>(paramList, factoryMapIn, factoryManagersIn);
145 if (factoryName ==
"BlockedRAPFactory")
return BuildRAPFactory<BlockedRAPFactory>(paramList, factoryMapIn, factoryManagersIn);
146 if (factoryName ==
"BrickAggregationFactory")
return Build2<BrickAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
147 if (factoryName ==
"ClassicalMapFactory")
return Build2<ClassicalMapFactory>(paramList, factoryMapIn, factoryManagersIn);
148 if (factoryName ==
"ClassicalPFactory")
return Build2<ClassicalPFactory>(paramList, factoryMapIn, factoryManagersIn);
149 if (factoryName ==
"CloneRepartitionInterface")
return Build2<CloneRepartitionInterface>(paramList, factoryMapIn, factoryManagersIn);
150 if (factoryName ==
"CoarseMapFactory")
return Build2<CoarseMapFactory>(paramList, factoryMapIn, factoryManagersIn);
151 if (factoryName ==
"CoarseningVisualizationFactory")
return Build2<CoarseningVisualizationFactory>(paramList, factoryMapIn, factoryManagersIn);
152 if (factoryName ==
"CoalesceDropFactory")
return Build2<CoalesceDropFactory>(paramList, factoryMapIn, factoryManagersIn);
153 if (factoryName ==
"SmooVecCoalesceDropFactory")
return Build2<SmooVecCoalesceDropFactory>(paramList, factoryMapIn, factoryManagersIn);
154 if (factoryName ==
"ConstraintFactory")
return Build2<ConstraintFactory>(paramList, factoryMapIn, factoryManagersIn);
155 if (factoryName ==
"CoordinatesTransferFactory")
return Build2<CoordinatesTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
156 if (factoryName ==
"DirectSolver")
return BuildDirectSolver(paramList, factoryMapIn, factoryManagersIn);
157 if (factoryName ==
"DropNegativeEntriesFactory")
return Build2<DropNegativeEntriesFactory>(paramList, factoryMapIn, factoryManagersIn);
158 if (factoryName ==
"EminPFactory")
return Build2<EminPFactory>(paramList, factoryMapIn, factoryManagersIn);
159 if (factoryName ==
"FilteredAFactory")
return Build2<FilteredAFactory>(paramList, factoryMapIn, factoryManagersIn);
160 if (factoryName ==
"FineLevelInputDataFactory")
return Build2<FineLevelInputDataFactory>(paramList, factoryMapIn, factoryManagersIn);
161 if (factoryName ==
"GeneralGeometricPFactory")
return Build2<GeneralGeometricPFactory>(paramList, factoryMapIn, factoryManagersIn);
162 if (factoryName ==
"ReplicatePFactory")
return Build2<ReplicatePFactory>(paramList, factoryMapIn, factoryManagersIn);
163 if (factoryName ==
"CombinePFactory")
return Build2<CombinePFactory>(paramList, factoryMapIn, factoryManagersIn);
164 if (factoryName ==
"GenericRFactory")
return Build2<GenericRFactory>(paramList, factoryMapIn, factoryManagersIn);
165 if (factoryName ==
"GeometricInterpolationPFactory")
return Build2<GeometricInterpolationPFactory>(paramList, factoryMapIn, factoryManagersIn);
166 if (factoryName ==
"HybridAggregationFactory")
return Build2<HybridAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
167 if (factoryName ==
"InterfaceAggregationFactory")
return Build2<InterfaceAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
168 if (factoryName ==
"InterfaceMappingTransferFactory")
return Build2<InterfaceMappingTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
169 if (factoryName ==
"InverseApproximationFactory")
return Build2<InverseApproximationFactory>(paramList, factoryMapIn, factoryManagersIn);
170 if (factoryName ==
"InitialBlockNumberFactory")
return Build2<InitialBlockNumberFactory>(paramList, factoryMapIn, factoryManagersIn);
171 if (factoryName ==
"LineDetectionFactory")
return Build2<LineDetectionFactory>(paramList, factoryMapIn, factoryManagersIn);
174 if (factoryName ==
"MapTransferFactory")
return Build2<MapTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
175 if (factoryName ==
"MatrixAnalysisFactory")
return Build2<MatrixAnalysisFactory>(paramList, factoryMapIn, factoryManagersIn);
176 if (factoryName ==
"MultiVectorTransferFactory")
return Build2<MultiVectorTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
179 if (factoryName ==
"NotayAggregationFactory")
return Build2<NotayAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
180 if (factoryName ==
"NullspaceFactory")
return Build2<NullspaceFactory>(paramList, factoryMapIn, factoryManagersIn);
181 if (factoryName ==
"NullspacePresmoothFactory")
return Build2<NullspacePresmoothFactory>(paramList, factoryMapIn, factoryManagersIn);
182 if (factoryName ==
"PatternFactory")
return Build2<PatternFactory>(paramList, factoryMapIn, factoryManagersIn);
183 if (factoryName ==
"PgPFactory")
return Build2<PgPFactory>(paramList, factoryMapIn, factoryManagersIn);
184 if (factoryName ==
"SaPFactory")
return Build2<SaPFactory>(paramList, factoryMapIn, factoryManagersIn);
185 if (factoryName ==
"RAPFactory")
return BuildRAPFactory<RAPFactory>(paramList, factoryMapIn, factoryManagersIn);
186 if (factoryName ==
"RAPShiftFactory")
return BuildRAPFactory<RAPShiftFactory>(paramList, factoryMapIn, factoryManagersIn);
187 if (factoryName ==
"RebalanceAcFactory")
return Build2<RebalanceAcFactory>(paramList, factoryMapIn, factoryManagersIn);
188 if (factoryName ==
"RebalanceTransferFactory")
return Build2<RebalanceTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
189 if (factoryName ==
"RegionRFactory")
return Build2<RegionRFactory>(paramList, factoryMapIn, factoryManagersIn);
190 if (factoryName ==
"RegionRFactory_kokkos")
return Build2<RegionRFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
191 if (factoryName ==
"ReorderBlockAFactory")
return Build2<ReorderBlockAFactory>(paramList, factoryMapIn, factoryManagersIn);
192 if (factoryName ==
"RepartitionInterface")
return Build2<RepartitionInterface>(paramList, factoryMapIn, factoryManagersIn);
193 if (factoryName ==
"ScaledNullspaceFactory")
return Build2<ScaledNullspaceFactory>(paramList, factoryMapIn, factoryManagersIn);
194 if (factoryName ==
"SegregatedAFactory")
return Build2<SegregatedAFactory>(paramList, factoryMapIn, factoryManagersIn);
195 if (factoryName ==
"SemiCoarsenPFactory")
return Build2<SemiCoarsenPFactory>(paramList, factoryMapIn, factoryManagersIn);
196 if (factoryName ==
"StructuredAggregationFactory")
return Build2<StructuredAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
197 if (factoryName ==
"StructuredLineDetectionFactory")
return Build2<StructuredLineDetectionFactory>(paramList, factoryMapIn, factoryManagersIn);
198 if (factoryName ==
"SubBlockAFactory")
return Build2<SubBlockAFactory>(paramList, factoryMapIn, factoryManagersIn);
199 if (factoryName ==
"TentativePFactory")
return Build2<TentativePFactory>(paramList, factoryMapIn, factoryManagersIn);
200 if (factoryName ==
"ToggleCoordinatesTransferFactory")
return BuildToggleCoordinatesTransferFactory(paramList, factoryMapIn, factoryManagersIn);
201 if (factoryName ==
"TogglePFactory")
return BuildTogglePFactory<TogglePFactory>(paramList, factoryMapIn, factoryManagersIn);
202 if (factoryName ==
"TransPFactory")
return Build2<TransPFactory>(paramList, factoryMapIn, factoryManagersIn);
203 if (factoryName ==
"RfromP_Or_TransP")
return Build2<RfromP_Or_TransP>(paramList, factoryMapIn, factoryManagersIn);
204 if (factoryName ==
"TrilinosSmoother")
return BuildTrilinosSmoother(paramList, factoryMapIn, factoryManagersIn);
205 if (factoryName ==
"UncoupledAggregationFactory")
return Build2<UncoupledAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
206 if (factoryName ==
"UnsmooshFactory")
return Build2<UnsmooshFactory>(paramList, factoryMapIn, factoryManagersIn);
207 if (factoryName ==
"UserAggregationFactory")
return Build2<UserAggregationFactory>(paramList, factoryMapIn, factoryManagersIn);
208 if (factoryName ==
"UserPFactory")
return Build2<UserPFactory>(paramList, factoryMapIn, factoryManagersIn);
209 if (factoryName ==
"VariableDofLaplacianFactory")
return Build2<VariableDofLaplacianFactory>(paramList, factoryMapIn, factoryManagersIn);
210 if (factoryName ==
"ZeroSubBlockAFactory")
return Build2<ZeroSubBlockAFactory>(paramList, factoryMapIn, factoryManagersIn);
211 if (factoryName ==
"CoalesceDropFactory_kokkos")
return Build2<CoalesceDropFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
212 if (factoryName ==
"GeometricInterpolationPFactory_kokkos")
return Build2<GeometricInterpolationPFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
213 #ifdef HAVE_MUELU_DEPRECATED_CODE
214 if (factoryName ==
"SaPFactory_kokkos")
return Build2<SaPFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
216 if (factoryName ==
"SemiCoarsenPFactory_kokkos")
return Build2<SemiCoarsenPFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
217 if (factoryName ==
"StructuredAggregationFactory_kokkos")
return Build2<StructuredAggregationFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
218 if (factoryName ==
"TentativePFactory_kokkos")
return Build2<TentativePFactory_kokkos>(paramList, factoryMapIn, factoryManagersIn);
219 if (factoryName ==
"MatrixFreeTentativePFactory")
return Build2<MatrixFreeTentativePFactory>(paramList, factoryMapIn, factoryManagersIn);
222 if (factoryName ==
"CoarseMapFactory_kokkos")
return Build2<CoarseMapFactory>(paramList, factoryMapIn, factoryManagersIn);
223 if (factoryName ==
"CoordinatesTransferFactory_kokkos")
return Build2<CoordinatesTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
225 if (factoryName ==
"ZoltanInterface") {
226 #if defined(HAVE_MUELU_ZOLTAN) && defined(HAVE_MPI)
227 return Build2<ZoltanInterface>(paramList, factoryMapIn, factoryManagersIn);
230 #endif // HAVE_MUELU_ZOLTAN && HAVE_MPI
232 if (factoryName ==
"Zoltan2Interface") {
233 #if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
234 return Build2<Zoltan2Interface>(paramList, factoryMapIn, factoryManagersIn);
237 #endif // HAVE_MUELU_ZOLTAN2 && HAVE_MPI
239 if (factoryName ==
"IsorropiaInterface") {
240 #if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
241 return Build2<IsorropiaInterface>(paramList, factoryMapIn, factoryManagersIn);
244 #endif // HAVE_MUELU_ZOLTAN2 && HAVE_MPI
247 if (factoryName ==
"NodePartitionInterface") {
248 #if defined(HAVE_MPI)
249 return Build2<NodePartitionInterface>(paramList, factoryMapIn, factoryManagersIn);
255 if (factoryName ==
"RepartitionFactory") {
257 return Build2<RepartitionFactory>(paramList, factoryMapIn, factoryManagersIn);
262 if (factoryName ==
"RepartitionHeuristicFactory") {
264 return Build2<RepartitionHeuristicFactory>(paramList, factoryMapIn, factoryManagersIn);
270 if (factoryName ==
"BlockedCoordinatesTransferFactory")
return BuildBlockedCoordFactory<BlockedCoordinatesTransferFactory>(paramList, factoryMapIn, factoryManagersIn);
271 if (factoryName ==
"BlockedDirectSolver")
return BuildBlockedDirectSolver(paramList, factoryMapIn, factoryManagersIn);
272 if (factoryName ==
"BlockedGaussSeidelSmoother")
return BuildBlockedSmoother<BlockedGaussSeidelSmoother>(paramList, factoryMapIn, factoryManagersIn);
273 if (factoryName ==
"BlockedJacobiSmoother")
return BuildBlockedSmoother<BlockedJacobiSmoother>(paramList, factoryMapIn, factoryManagersIn);
274 if (factoryName ==
"BlockedPFactory")
return BuildBlockedFactory<BlockedPFactory>(paramList, factoryMapIn, factoryManagersIn);
275 if (factoryName ==
"BraessSarazinSmoother")
return BuildBlockedSmoother<BraessSarazinSmoother>(paramList, factoryMapIn, factoryManagersIn);
276 if (factoryName ==
"IndefiniteBlockDiagonalSmoother")
return BuildBlockedSmoother<IndefBlockedDiagonalSmoother>(paramList, factoryMapIn, factoryManagersIn);
277 if (factoryName ==
"SimpleSmoother")
return BuildBlockedSmoother<SimpleSmoother>(paramList, factoryMapIn, factoryManagersIn);
278 if (factoryName ==
"SchurComplementFactory")
return Build2<SchurComplementFactory>(paramList, factoryMapIn, factoryManagersIn);
279 if (factoryName ==
"RebalanceBlockRestrictionFactory")
return BuildBlockedFactory<RebalanceBlockRestrictionFactory>(paramList, factoryMapIn, factoryManagersIn);
280 if (factoryName ==
"RebalanceBlockAcFactory")
return BuildBlockedFactory<RebalanceBlockAcFactory>(paramList, factoryMapIn, factoryManagersIn);
281 if (factoryName ==
"RebalanceBlockInterpolationFactory")
return BuildBlockedFactory<RebalanceBlockInterpolationFactory>(paramList, factoryMapIn, factoryManagersIn);
283 if (factoryName ==
"RepartitionBlockDiagonalFactory")
return Build2<RepartitionBlockDiagonalFactory>(paramList, factoryMapIn, factoryManagersIn);
285 #ifdef HAVE_MUELU_TEKO
286 if (factoryName ==
"TekoSmoother")
return BuildTekoSmoother(paramList, factoryMapIn, factoryManagersIn);
288 if (factoryName ==
"UzawaSmoother")
return BuildBlockedSmoother<UzawaSmoother>(paramList, factoryMapIn, factoryManagersIn);
291 #ifdef HAVE_MUELU_MATLAB
292 if (factoryName ==
"SingleLevelMatlabFactory")
return Build2<SingleLevelMatlabFactory>(paramList, factoryMapIn, factoryManagersIn);
293 if (factoryName ==
"TwoLevelMatlabFactory")
return Build2<TwoLevelMatlabFactory>(paramList, factoryMapIn, factoryManagersIn);
294 if (factoryName ==
"MatlabSmoother")
return BuildMatlabSmoother(paramList, factoryMapIn, factoryManagersIn);
297 #ifdef HAVE_MUELU_INTREPID2
298 if (factoryName ==
"IntrepidPCoarsenFactory")
return Build2<IntrepidPCoarsenFactory>(paramList, factoryMapIn, factoryManagersIn);
302 if (factoryMapIn.find(factoryName) != factoryMapIn.end()) {
304 "MueLu::FactoryFactory: Error during the parsing of: " << std::endl
305 << paramList << std::endl
306 <<
"'" << factoryName <<
"' is not a factory name but an existing instance of a factory." << std::endl
307 <<
"Extra parameters cannot be specified after the creation of the object." << std::endl
309 <<
"Correct syntaxes includes:" << std::endl
310 <<
" <Parameter name=\"...\" type=\"string\" value=\"" << factoryName <<
"\"/>" << std::endl
312 <<
" <ParameterList name=\"...\"><Parameter name=\"factory\" type=\"string\" value=\"" << factoryName <<
"\"/></ParameterList>" << std::endl);
314 return factoryMapIn.find(factoryName)->second;
322 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
327 const char* strarray[] = {
"A",
"P",
"R",
"Graph",
"UnAmalgamationInfo",
"Aggregates",
"Nullspace",
"TransferFactory",
"DofsPerNode"};
328 #define arraysize(ar) (sizeof(ar) / sizeof(ar[0]))
329 std::vector<std::string> v(strarray, strarray +
arraysize(strarray));
330 for (
size_t i = 0; i < v.size(); ++i)
332 factory->SetFactory(v[i], BuildFactory(paramList.
getEntry(v[i]), factoryMapIn, factoryManagersIn));
337 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
348 const std::string& pName = validParamList->
name(param);
358 paramListWithFactories.
set(pName, generatingFact);
360 if (pName ==
"ParameterList") {
366 paramListWithFactories.
set(pName, subList);
374 factory->SetParameterList(paramListWithFactories);
379 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
383 if (paramList.
isSublist(
"TransferFactories") ==
false) {
384 factory = Build2<T>(paramList, factoryMapIn, factoryManagersIn);
390 paramListNonConst->
remove(
"TransferFactories");
392 factory = Build2<T>(*paramListNonConst, factoryMapIn, factoryManagersIn);
396 factory->AddTransferFactory(p);
403 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
407 if (paramList.
isSublist(
"TransferFactories") ==
false) {
409 factory = Build2<T>(paramList, factoryMapIn, factoryManagersIn);
415 paramListNonConst->
remove(
"TransferFactories");
418 factory = Build2<T>(*paramListNonConst, factoryMapIn, factoryManagersIn);
422 int numProlongatorFactories = 0;
423 int numPtentFactories = 0;
424 int numCoarseNspFactories = 0;
426 size_t foundNsp = transferFactories->
name(param).find(
"Nullspace");
427 if (foundNsp != std::string::npos && foundNsp == 0 && transferFactories->
name(param).length() == 10) {
428 numCoarseNspFactories++;
431 size_t foundPtent = transferFactories->
name(param).find(
"Ptent");
432 if (foundPtent != std::string::npos && foundPtent == 0 && transferFactories->
name(param).length() == 6) {
436 size_t foundP = transferFactories->
name(param).find(
"P");
437 if (foundP != std::string::npos && foundP == 0 && transferFactories->
name(param).length() == 2) {
438 numProlongatorFactories++;
444 TEUCHOS_TEST_FOR_EXCEPTION(numProlongatorFactories < 2,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleP: The TogglePFactory needs at least two different prolongation operators. The factories have to be provided using the names P%i and Nullspace %i, where %i denotes a number between 1 and 9.");
447 std::vector<Teuchos::ParameterEntry> prolongatorFactoryNames(numProlongatorFactories);
448 std::vector<Teuchos::ParameterEntry> coarseNspFactoryNames(numProlongatorFactories);
449 std::vector<Teuchos::ParameterEntry> ptentFactoryNames(numProlongatorFactories);
452 size_t foundNsp = transferFactories->
name(param).find(
"Nullspace");
453 if (foundNsp != std::string::npos && foundNsp == 0 && transferFactories->
name(param).length() == 10) {
454 int number = atoi(&(transferFactories->
name(param).at(9)));
455 TEUCHOS_TEST_FOR_EXCEPTION(number < 1 || number > numProlongatorFactories,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleP: Please use the format Nullspace%i with %i an integer between 1 and the maximum number of prolongation operators in TogglePFactory!");
456 coarseNspFactoryNames[number - 1] = transferFactories->
entry(param);
459 size_t foundPtent = transferFactories->
name(param).find(
"Ptent");
460 if (foundPtent != std::string::npos && foundPtent == 0 && transferFactories->
name(param).length() == 6) {
461 int number = atoi(&(transferFactories->
name(param).at(5)));
462 TEUCHOS_TEST_FOR_EXCEPTION(number < 1 || number > numPtentFactories,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleP: Please use the format Ptent%i with %i an integer between 1 and the maximum number of prolongation operators in TogglePFactory!");
463 ptentFactoryNames[number - 1] = transferFactories->
entry(param);
466 size_t foundP = transferFactories->
name(param).find(
"P");
467 if (foundP != std::string::npos && foundP == 0 && transferFactories->
name(param).length() == 2) {
468 int number = atoi(&(transferFactories->
name(param).at(1)));
469 TEUCHOS_TEST_FOR_EXCEPTION(number < 1 || number > numProlongatorFactories,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleP: Please use the format P%i with %i an integer between 1 and the maximum number of prolongation operators in TogglePFactory!");
470 prolongatorFactoryNames[number - 1] = transferFactories->
entry(param);
476 for (std::vector<Teuchos::ParameterEntry>::const_iterator it = prolongatorFactoryNames.begin(); it != prolongatorFactoryNames.end(); ++it) {
478 factory->AddProlongatorFactory(p);
482 for (std::vector<Teuchos::ParameterEntry>::const_iterator it = ptentFactoryNames.begin(); it != ptentFactoryNames.end(); ++it) {
484 factory->AddPtentFactory(p);
488 for (std::vector<Teuchos::ParameterEntry>::const_iterator it = coarseNspFactoryNames.begin(); it != coarseNspFactoryNames.end(); ++it) {
490 factory->AddCoarseNullspaceFactory(p);
496 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
499 TEUCHOS_TEST_FOR_EXCEPTION(paramList.
isSublist(
"TransferFactories") ==
false,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleCoordinatesTransferFactory: the ToggleCoordinatesTransferFactory needs a sublist 'TransferFactories' containing information about the subfactories for coordinate transfer!");
503 paramListNonConst->remove(
"TransferFactories");
506 factory = Build2<ToggleCoordinatesTransferFactory>(*paramListNonConst, factoryMapIn, factoryManagersIn);
510 int numCoordTransferFactories = 0;
512 size_t foundCoordinates = transferFactories->
name(param).find(
"Coordinates");
513 if (foundCoordinates != std::string::npos && foundCoordinates == 0 && transferFactories->
name(param).length() == 12) {
514 numCoordTransferFactories++;
518 TEUCHOS_TEST_FOR_EXCEPTION(numCoordTransferFactories != 2,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleCoordinatesTransfer: The ToggleCoordinatesTransferFactory needs two (different) coordinate transfer factories. The factories have to be provided using the names Coordinates%i, where %i denotes a number between 1 and 9.");
521 std::vector<Teuchos::ParameterEntry> coarseCoordsFactoryNames(numCoordTransferFactories);
524 size_t foundCoords = transferFactories->
name(param).find(
"Coordinates");
525 if (foundCoords != std::string::npos && foundCoords == 0 && transferFactories->
name(param).length() == 12) {
526 int number = atoi(&(transferFactories->
name(param).at(11)));
527 TEUCHOS_TEST_FOR_EXCEPTION(number < 1 || number > numCoordTransferFactories,
Exceptions::RuntimeError,
"FactoryFactory::BuildToggleCoordinatesTransfer: Please use the format Coordinates%i with %i an integer between 1 and the maximum number of coordinate transfer factories in ToggleCoordinatesTransferFactory!");
528 coarseCoordsFactoryNames[number - 1] = transferFactories->
entry(param);
534 for (std::vector<Teuchos::ParameterEntry>::const_iterator it = coarseCoordsFactoryNames.begin(); it != coarseCoordsFactoryNames.end(); ++it) {
542 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
544 if (paramList.
begin() == paramList.
end())
552 std::string type =
"";
553 if (paramList.
isParameter(
"type")) type = paramList.
get<std::string>(
"type");
555 if (paramList.
isParameter(
"overlap")) overlap = paramList.
get<
int>(
"overlap");
569 if (paramList.
isParameter(
"LineDetection_Layers")) {
571 trilSmoo->
SetFactory(
"LineDetection_Layers", generatingFact);
573 if (paramList.
isParameter(
"LineDetection_VertLineIds")) {
575 trilSmoo->
SetFactory(
"LineDetection_Layers", generatingFact);
579 trilSmoo->
SetFactory(
"CoarseNumZLayers", generatingFact);
590 #ifdef HAVE_MUELU_MATLAB
591 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
594 if (paramList.
begin() == paramList.
end())
609 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
611 if (paramList.
begin() == paramList.
end())
617 if (paramList.
isParameter(
"type")) type = paramList.
get<std::string>(
"type");
625 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
632 std::vector<RCP<FactoryManager> > facManagers;
636 bool blockExists =
true;
637 while (blockExists ==
true) {
638 std::stringstream ss;
639 ss <<
"block" << blockid;
641 if (paramList.
isSublist(ss.str()) ==
true) {
649 std::string facManagerName = b->
get<std::string>(
"group");
659 M->SetFactory(b->
name(param), p);
664 M->SetIgnoreUserData(
true);
665 facManagers.push_back(M);
666 paramListNonConst->
remove(ss.str());
675 RCP<T> bs = Build2<T>(*paramListNonConst, factoryMapIn, factoryManagersIn);
683 for (
int i = 0; i < Teuchos::as<int>(facManagers.size()); i++) {
684 bs->AddFactoryManager(facManagers[i], i);
690 #ifdef HAVE_MUELU_TEKO
691 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
696 paramListNonConst->
remove(
"Inverse Factory Library");
699 RCP<TekoSmoother> bs = Build2<TekoSmoother>(*paramListNonConst, factoryMapIn, factoryManagersIn);
714 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
722 if (paramList.
isParameter(
"type")) type = paramList.
get<std::string>(
"type");
730 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
733 RCP<T> pfac = Teuchos::null;
739 std::vector<RCP<FactoryManager> > facManagers;
743 bool blockExists =
true;
744 while (blockExists ==
true) {
745 std::stringstream ss;
746 ss <<
"block" << blockid;
748 if (paramList.
isSublist(ss.str()) ==
true) {
756 std::string facManagerName = b->
get<std::string>(
"group");
766 M->SetFactory(b->
name(param), p);
771 M->SetIgnoreUserData(
true);
772 facManagers.push_back(M);
773 paramListNonConst->
remove(ss.str());
782 pfac = Build2<T>(*paramListNonConst, factoryMapIn, factoryManagersIn);
785 for (
size_t i = 0; i < facManagers.size(); i++) {
786 pfac->AddFactoryManager(facManagers[i]);
792 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
795 RCP<T> pfac = Teuchos::null;
801 std::vector<RCP<const FactoryBase> > facBase;
805 bool blockExists =
true;
806 while (blockExists ==
true) {
807 std::stringstream ss;
808 ss <<
"block" << blockid;
810 if (paramList.
isSublist(ss.str()) ==
true) {
817 facBase.push_back(p);
821 paramListNonConst->
remove(ss.str());
830 pfac = Build2<T>(*paramListNonConst, factoryMapIn, factoryManagersIn);
833 for (
size_t i = 0; i < facBase.size(); i++) {
834 pfac->AddFactory(facBase[i]);
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
RCP< ToggleCoordinatesTransferFactory > BuildToggleCoordinatesTransferFactory(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
const std::string & name() const
This class specifies the default factory that should generate some data on a Level if the data does n...
ParameterList & setEntry(const std::string &name, U &&entry)
ConstIterator end() const
virtual RCP< const FactoryBase > BuildFactory(const Teuchos::ParameterEntry ¶m, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
: Interpret Factory parameter list and build new factory
T & get(const std::string &name, T def_value)
Class that encapsulates external library smoothers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< T > Build2(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
RCP< T > BuildTogglePFactory(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
Ordinal numParams() const
RCP< T > BuildBlockedFactory(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void AddCoordTransferFactory(const RCP< const FactoryBase > &factory)
Add a coordinate transfer factory in the end of list of coordinate transfer factories.
RCP< FactoryBase > BuildDirectSolver(const Teuchos::ParameterList ¶mList, const FactoryMap &, const FactoryManagerMap &) const
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
void SetTekoParameters(RCP< ParameterList > tekoParams)
bool isParameter(const std::string &name) const
std::map< std::string, RCP< FactoryManagerBase > > FactoryManagerMap
bool remove(std::string const &name, bool throwIfNotExists=true)
RCP< FactoryBase > BuildTrilinosSmoother(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
TrilinosSmoother.
RCP< T > Build(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
virtual void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
void SetSmootherPrototypes(RCP< SmootherPrototype > preAndPostSmootherPrototype)
Set smoother prototypes.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
std::map< std::string, RCP< const FactoryBase > > FactoryMap
bool isSublist(const std::string &name) const
RCP< FactoryBase > BuildBlockedDirectSolver(const Teuchos::ParameterList ¶mList, const FactoryMap &, const FactoryManagerMap &) const
params_t::ConstIterator ConstIterator
static const RCP< const NoFactory > getRCP()
Static Get() functions.
ConstIterator begin() const
direct solver for nxn blocked matrices
const ParameterEntry & entry(ConstIterator i) const
Class that encapsulates Matlab smoothers.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
RCP< FactoryBase > BuildBlockedSmoother(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
RCP< T > BuildRAPFactory(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
bool isType(const std::string &name) const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Exception throws to report errors in the internal logical of the program.
RCP< FactoryBase > BuildTekoSmoother(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
ParameterEntry & getEntry(const std::string &name)
RCP< FactoryBase > BuildMatlabSmoother(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
MatlabSmoother.
RCP< T > BuildBlockedCoordFactory(const Teuchos::ParameterList ¶mList, const FactoryMap &factoryMapIn, const FactoryManagerMap &factoryManagersIn) const
void SetFactory(const std::string &varName, const RCP< const FactoryBase > &factory)
Custom SetFactory.