46 #ifndef MUELU_TWOLEVELMATLABFACTORY_DEF_HPP
47 #define MUELU_TWOLEVELMATLABFACTORY_DEF_HPP
48 #include <Xpetra_Matrix.hpp>
49 #include <Xpetra_MultiVector.hpp>
52 #include "MueLu_Aggregates.hpp"
53 #include "MueLu_AmalgamationInfo.hpp"
60 #ifdef HAVE_MUELU_MATLAB
65 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
67 : hasDeclaredInput_(false) { }
69 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
72 validParamList->
set<std::string>(
"Provides" ,
"" ,
"A comma-separated list of objects provided on the coarse level by the TwoLevelMatlabFactory");
73 validParamList->
set<std::string>(
"Needs Fine" ,
"",
"A comma-separated list of objects needed on the fine level by the TwoLevelMatlabFactory");
74 validParamList->
set<std::string>(
"Needs Coarse" ,
"",
"A comma-separated list of objects needed on the coarse level by the TwoLevelMatlabFactory");
75 validParamList->
set<std::string>(
"Function" ,
"" ,
"The name of the Matlab MEX function to call for Build()");
76 return validParamList;
79 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
84 const std::string str_nf = pL.get<std::string>(
"Needs Fine");
85 const std::string str_nc = pL.get<std::string>(
"Needs Coarse");
88 for(
auto fineNeed : needsFine_)
91 this->Input(fineLevel, fineNeed);
93 for(
auto coarseNeed : needsCoarse_)
96 this->Input(coarseLevel, coarseNeed);
98 hasDeclaredInput_ =
true;
101 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
109 string needsFine = pL.
get<
string>(
"Needs Fine");
110 string needsCoarse = pL.
get<
string>(
"Needs Coarse");
111 vector<RCP<MuemexArg>> InputArgs = processNeeds<Scalar, LocalOrdinal, GlobalOrdinal, Node>(
this, needsFine, fineLevel);
112 vector<RCP<MuemexArg>> InputArgsCoarse = processNeeds<Scalar, LocalOrdinal, GlobalOrdinal, Node>(
this, needsCoarse, coarseLevel);
114 InputArgs.reserve(InputArgs.size() + InputArgsCoarse.size());
115 InputArgs.insert(InputArgs.begin(), InputArgsCoarse.begin(), InputArgsCoarse.end());
118 string provides = pL.
get<
string>(
"Provides");
121 string matlabFunction = pL.
get<
string>(
"Function");
122 if(!matlabFunction.length())
123 throw runtime_error(
"Invalid matlab function name");
124 vector<RCP<MuemexArg>> mexOutput =
callMatlab(matlabFunction, numProvides, InputArgs);
125 processProvides<Scalar, LocalOrdinal, GlobalOrdinal, Node>(mexOutput,
this, provides, coarseLevel);
128 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
130 std::ostringstream out;
132 out <<
"TwoLevelMatlabFactory["<<pL.
get<std::string>(
"Function")<<
"]";
139 #define MUELU_TWOLEVELMATLABFACTORY_SHORT
140 #endif // HAVE_MUELU_MATLAB
142 #endif // MUELU_TWOLEVELMATLABFACTORY_DEF_HPP
std::vector< std::string > tokenizeList(const std::string ¶ms)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
bool IsParamMuemexVariable(const std::string &name)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
Teuchos::RCP< Teuchos::ParameterList > getInputParamList()
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
std::string description() const
@ name Description
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
std::vector< RCP< MuemexArg > > callMatlab(std::string function, int numOutputs, std::vector< RCP< MuemexArg >> args)