MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MatlabUtils_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_MATLABUTILS_DECL_HPP
11#define MUELU_MATLABUTILS_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14
15#if !defined(HAVE_MUELU_MATLAB) || !defined(HAVE_MUELU_TPETRA)
16#error "Muemex requires MATLAB, Epetra and Tpetra."
17#else
18
19// Matlab fwd style declarations
20struct mxArray_tag;
21typedef struct mxArray_tag mxArray;
22typedef size_t mwIndex;
23
24#include <string>
25#include <complex>
26#include <stdexcept>
27#include <Teuchos_ParameterList.hpp>
28#include <Teuchos_RCP.hpp>
29#include <Teuchos_DefaultComm.hpp>
30#include "MueLu_Factory.hpp"
35#include "MueLu_Graph_fwd.hpp"
36#ifdef HAVE_MUELU_EPETRA
37#include "Epetra_MultiVector.h"
38#include "Epetra_CrsMatrix.h"
39#endif
40#include "Tpetra_CrsMatrix_decl.hpp"
41#ifdef HAVE_MUELU_EPETRA
42#include "Xpetra_EpetraCrsMatrix.hpp"
43#endif
44#include "Xpetra_MapFactory.hpp"
45#include "Xpetra_CrsGraph.hpp"
46#include "Xpetra_VectorFactory.hpp"
47#include <Tpetra_Core.hpp>
48
49#include "Kokkos_DynRankView.hpp"
50
51namespace MueLu {
52
81
82typedef Tpetra::KokkosCompat::KokkosDeviceWrapperNode<Kokkos::Serial, Kokkos::HostSpace> mm_node_t;
83typedef typename Tpetra::Map<>::local_ordinal_type mm_LocalOrd; // these are used for LocalOrdinal and GlobalOrdinal of all xpetra/tpetra templated types
84typedef typename Tpetra::Map<>::global_ordinal_type mm_GlobalOrd;
85typedef std::complex<double> complex_t;
86typedef Tpetra::Map<> muemex_map_type;
87typedef Tpetra::CrsMatrix<double, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Tpetra_CrsMatrix_double;
88typedef Tpetra::CrsMatrix<complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Tpetra_CrsMatrix_complex;
89typedef Tpetra::MultiVector<double, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Tpetra_MultiVector_double;
90typedef Tpetra::MultiVector<complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Tpetra_MultiVector_complex;
91typedef Xpetra::Map<mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_map;
92typedef Xpetra::Vector<mm_LocalOrd, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_ordinal_vector;
93typedef Xpetra::Matrix<double, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_Matrix_double;
94typedef Xpetra::Matrix<complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_Matrix_complex;
95typedef Xpetra::CrsGraph<mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_CrsGraph;
96typedef Xpetra::MultiVector<double, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_MultiVector_double;
97typedef Xpetra::MultiVector<complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t> Xpetra_MultiVector_complex;
103
104#ifdef HAVE_MUELU_INTREPID2
105typedef Kokkos::DynRankView<mm_LocalOrd, typename mm_node_t::device_type> FieldContainer_ordinal;
106#endif
107
109 public:
110 MuemexArg(MuemexType dataType) { type = dataType; }
112};
113
114template <typename T>
115MuemexType getMuemexType(const T& data);
116
117template <typename T>
118class MuemexData : public MuemexArg {
119 public:
120 MuemexData(T& data); // Construct from pre-existing data, to pass to MATLAB.
121 MuemexData(T& data, MuemexType type); // Construct from pre-existing data, to pass to MATLAB.
122 MuemexData(const mxArray* mxa); // Construct from MATLAB array, to get from MATLAB.
123 mxArray* convertToMatlab(); // Create a MATLAB object and copy this data to it
124 T& getData(); // Set and get methods
125 void setData(T& data);
126
127 private:
129};
130
131template <typename T>
132MuemexType getMuemexType(const T& data);
133
134template <typename T>
136
137template <typename T>
139
140template <typename T>
142
143// Add data to level. Set the keep flag on the data to "user-provided" so it's not deleted.
144template <typename T>
145void addLevelVariable(const T& data, std::string& name, Level& lvl, const FactoryBase* fact = NoFactory::get());
146
147template <typename T>
148const T& getLevelVariable(std::string& name, Level& lvl);
149
150// Functions used to put data through matlab factories - first arg is "this" pointer of matlab factory
151template <typename Scalar = double, typename LocalOrdinal = mm_LocalOrd, typename GlobalOrdinal = mm_GlobalOrd, typename Node = mm_node_t>
152std::vector<Teuchos::RCP<MuemexArg> > processNeeds(const Factory* factory, std::string& needsParam, Level& lvl);
153
154template <typename Scalar = double, typename LocalOrdinal = mm_LocalOrd, typename GlobalOrdinal = mm_GlobalOrd, typename Node = mm_node_t>
155void processProvides(std::vector<Teuchos::RCP<MuemexArg> >& mexOutput, const Factory* factory, std::string& providesParam, Level& lvl);
156
157// create a sparse array in Matlab
158template <typename Scalar>
159mxArray* createMatlabSparse(int numRows, int numCols, int nnz);
160template <typename Scalar>
161mxArray* createMatlabMultiVector(int numRows, int numCols);
162template <typename Scalar>
163void fillMatlabArray(Scalar* array, const mxArray* mxa, int n);
164int* mwIndex_to_int(int N, mwIndex* mwi_array);
165bool isValidMatlabAggregates(const mxArray* mxa);
166bool isValidMatlabGraph(const mxArray* mxa);
167std::vector<std::string> tokenizeList(const std::string& param);
168// The two callback functions that MueLu can call to run anything in MATLAB
169void callMatlabNoArgs(std::string function);
170std::vector<Teuchos::RCP<MuemexArg> > callMatlab(std::string function, int numOutputs, std::vector<Teuchos::RCP<MuemexArg> > args);
171Teuchos::RCP<Teuchos::ParameterList> getInputParamList();
172Teuchos::RCP<MuemexArg> convertMatlabVar(const mxArray* mxa);
173
174// trim from start
175static inline std::string& ltrim(std::string& s) {
176 s.erase(0, s.find_first_not_of(" "));
177 return s;
178}
179
180// trim from end
181static inline std::string& rtrim(std::string& s) {
182 s.erase(s.find_last_not_of(" "), std::string::npos);
183 return s;
184}
185
186// trim from both ends
187static inline std::string& trim(std::string& s) {
188 return ltrim(rtrim(s));
189}
190
191} // namespace MueLu
192
193#endif // HAVE_MUELU_MATLAB error handler
194#endif // MUELU_MATLABUTILS_DECL_HPP guard
int mwIndex
struct mxArray_tag mxArray
size_t mwIndex
MueLu::DefaultScalar Scalar
Container class for aggregation information.
minimal container class for storing amalgamation information
Base class for factories (e.g., R, P, and A_coarse).
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Lightweight MueLu representation of a compressed row storage graph.
Class that holds all level-specific information.
MuemexArg(MuemexType dataType)
static const NoFactory * get()
Namespace for MueLu classes and methods.
static std::string & trim(std::string &s)
MuemexType getMuemexType()
Tpetra::CrsMatrix< double, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Tpetra_CrsMatrix_double
Xpetra::Matrix< complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_Matrix_complex
bool isValidMatlabGraph(const mxArray *mxa)
Tpetra::MultiVector< double, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Tpetra_MultiVector_double
Xpetra::Vector< mm_LocalOrd, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_ordinal_vector
Teuchos::RCP< Teuchos::ParameterList > getInputParamList()
MueLu::Hierarchy< complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Hierarchy_complex
mxArray * createMatlabMultiVector(int numRows, int numCols)
Teuchos::RCP< MuemexArg > convertMatlabVar(const mxArray *mxa)
Tpetra::Map ::local_ordinal_type mm_LocalOrd
Xpetra::MultiVector< double, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_MultiVector_double
Tpetra::KokkosCompat::KokkosDeviceWrapperNode< Kokkos::Serial, Kokkos::HostSpace > mm_node_t
Tpetra::MultiVector< complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Tpetra_MultiVector_complex
bool isValidMatlabAggregates(const mxArray *mxa)
std::vector< RCP< MuemexArg > > callMatlab(std::string function, int numOutputs, std::vector< RCP< MuemexArg > > args)
Xpetra::CrsGraph< mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_CrsGraph
int * mwIndex_to_int(int N, mwIndex *mwi_array)
std::vector< Teuchos::RCP< MuemexArg > > processNeeds(const Factory *factory, std::string &needsParam, Level &lvl)
void processProvides(std::vector< Teuchos::RCP< MuemexArg > > &mexOutput, const Factory *factory, std::string &providesParam, Level &lvl)
Xpetra::Matrix< double, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_Matrix_double
MueLu::Aggregates< mm_LocalOrd, mm_GlobalOrd, mm_node_t > MAggregates
Tpetra::Map muemex_map_type
Tpetra::Map ::global_ordinal_type mm_GlobalOrd
const T & getLevelVariable(std::string &name, Level &lvl)
Xpetra::Map< mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_map
void fillMatlabArray(Scalar *array, const mxArray *mxa, int n)
std::complex< double > complex_t
void addLevelVariable(const T &data, std::string &name, Level &lvl, const FactoryBase *fact=NoFactory::get())
MueLu::AmalgamationInfo< mm_LocalOrd, mm_GlobalOrd, mm_node_t > MAmalInfo
Tpetra::CrsMatrix< complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Tpetra_CrsMatrix_complex
std::vector< std::string > tokenizeList(const std::string &params)
Xpetra::MultiVector< complex_t, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Xpetra_MultiVector_complex
T loadDataFromMatlab(const mxArray *mxa)
void callMatlabNoArgs(std::string function)
mxArray * createMatlabSparse(int numRows, int numCols, int nnz)
static std::string & rtrim(std::string &s)
MueLu::LWGraph< mm_LocalOrd, mm_GlobalOrd, mm_node_t > MGraph
static std::string & ltrim(std::string &s)
template mxArray * saveDataToMatlab(bool &data)
MueLu::Hierarchy< double, mm_LocalOrd, mm_GlobalOrd, mm_node_t > Hierarchy_double