EpetraExt  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
EpetraExt::EpetraExt_MatlabEngine Class Reference

A class which provides data and command access to Matlab from Epetra. More...

#include <EpetraExt_MatlabEngine.h>

 EpetraExt_MatlabEngine (const Epetra_Comm &Comm)
 EpetraExt_MatlabEngine constructor which creates a MatlabEngine object with a connection to an instance of the application Matlab by starting a new Matlab process. More...
 
 ~EpetraExt_MatlabEngine ()
 EpetraExt_MatlabEngine destructor which closes the connection to Matlab which causes the Matlab process to also exit. More...
 
int EvalString (char *command, char *outputBuffer=NULL, int outputBufferSize=-1)
 Sends a command to Matlab. More...
 
int PutMultiVector (const Epetra_MultiVector &A, const char *variableName)
 Puts a copy of the serial or distributed MultiVector into the Matlab workspace. More...
 
int PutRowMatrix (const Epetra_RowMatrix &A, const char *variableName, bool transA)
 Puts a copy of the serial or distributed RowMatrix into the Matlab workspace. More...
 
int PutCrsGraph (const Epetra_CrsGraph &A, const char *variableName, bool transA)
 not implemented yet More...
 
int PutSerialDenseMatrix (const Epetra_SerialDenseMatrix &A, const char *variableName, int proc=0)
 Puts a copy of the SerialDenseMatrix into the Matlab workspace. More...
 
int PutIntSerialDenseMatrix (const Epetra_IntSerialDenseMatrix &A, const char *variableName, int proc=0)
 Puts a copy of the IntSerialDenseMatrix into the Matlab workspace. More...
 
int PutBlockMap (const Epetra_BlockMap &blockMap, const char *variableName, bool transA)
 Puts a copy of the BlockMap or Map into the Matlab workspace. More...
 
int PutIntoMatlab (const char *variableName, mxArray *matlabA)
 Puts a mxArray into Matlab. More...
 
int GetMultiVector (const char *variableName, Epetra_MultiVector &A)
 Puts a Matlab variable into a MultiVector. More...
 
int GetSerialDenseMatrix (const char *variableName, Epetra_SerialDenseMatrix &A, int proc=0)
 Puts a Matlab variable into a SerialDenseMatrix on the specified PE. More...
 
int GetIntSerialDenseMatrix (const char *variableName, Epetra_IntSerialDenseMatrix &A, int proc=0)
 Puts a Matlab variable into a IntSerialDenseMatrix on the specified PE. More...
 
int GetCrsMatrix (const char *variableName, Epetra_CrsMatrix &A, bool getTrans)
 Puts a Matlab variable into a CrsMatrix. More...
 
int GetmxArrayDimensions (mxArray *matlabA, bool &isSparse, int &numRows, int &numCols, int &numNonZeros)
 Get general information about the mxArray. For internal use but can be used by an advanced user. More...
 
int GetmxArray (const char *variableName, mxArray **matlabA)
 Get a mxArray from Matlab. For internal use but can be used by an advanced user. More...
 

Detailed Description

A class which provides data and command access to Matlab from Epetra.

The EpetraExt_MatlabEngine class allows Epetra data objects to be

exported to Matlab and then operated on within Matlab using Matlab commands.

When an EpetraExt_MatlabEngine object is constructed a new instance of the application Matlab is started for EpetraExt_MatlabEngine to communicate with. All communication between EpetraExt_MatlabEngine and Matlab occurs on the root node (0) only. For parallel environments all objects are collected onto the root node before being put into Matlab. Object data is put into a mxArray which is then sent to the Matlab process. All objects passed to Matlab are copied into the Matlab memory space. So at the point when Matlab receives its copy of the data there is two copies of the mxArray in memory, one in the Epetra application and one in Matlab. Since Matlab has its own memory space the mxArray in the Epetra application should be deleted in order to free up memory space. All methods in EpetraExt_MatlabEngine that put Epetra objects into Matlab delete the temporary mxArray as soon as it is put into Matlab. If a user desires to create his/her own mxArray's and then send them to Matlab the method PutIntoMatlab can be used. It is important to note that PutIntoMatlab does NOT delete the mxArray it is passed. When the EpetraExt_MatlabEngine deconstructor is called the instance of Matlab that was started during construction of the EpetraExt_MatlabEngine object exits.

Error Codes

Build Instructions
These instructions can be found in the file Trilinos/packages/epetraext/doc/matlab.README.

Example Code
The following example code generates simple Epetra objects and then puts them into Matlab.

The point of this example is to illustrate the flow of calls when using EpetraExt_MatlabEngine. This example program can be found in the file Trilinos/packages/epetraext/example/matlab/cxx_main.cpp.

/*
//@HEADER
// ***********************************************************************
//
// AztecOO: An Object-Oriented Aztec Linear Solver Package
// Copyright (2002) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ***********************************************************************
//@HEADER
*/
#include "AztecOO_config.h"
#ifdef HAVE_MPI
#include "mpi.h"
#include "Epetra_MpiComm.h"
#else
#include "Epetra_SerialComm.h"
#endif
#include "Epetra_Map.h"
#include "Epetra_Vector.h"
#include "Epetra_CrsMatrix.h"
#include "AztecOO_Operator.h"
#include "Epetra_InvOperator.h"
#include <string>
// prototypes
int checkValues( double x, double y, string message = "", bool verbose = false) {
if (fabs((x-y)/x) > 0.01) {
return(1);
if (verbose) cout << "********** " << message << " check failed.********** " << endl;
}
else {
if (verbose) cout << message << " check OK." << endl;
return(0);
}
}
int main(int argc, char *argv[]) {
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
Epetra_MpiComm comm (MPI_COMM_WORLD);
#else
#endif
int MyPID = comm.MyPID();
bool verbose = false;
bool verbose1 = false;
// Check if we should print results to standard out
if (argc > 1) {
if ((argv[1][0] == '-') && (argv[1][1] == 'v')) {
verbose1 = true;
if (MyPID==0) verbose = true;
}
}
if (verbose1) cout << comm << endl;
// Uncomment the next three lines to debug in mpi mode
//int tmp;
//if (MyPID==0) cin >> tmp;
//comm.Barrier();
EPETRA_CHK_ERR(EpetraExt::MatlabFileToCrsMatrix("A.dat", comm, A));
x.Random();
A->Apply(x,b); // Generate RHS from x
Epetra_Vector xx(x); // Copy x to xx for later use
Epetra_LinearProblem problem(A, &x, &b);
// Construct a solver object for this problem
AztecOO solver(problem);
solver.SetAztecOption(AZ_precond, AZ_none);
if (!verbose1) solver.SetAztecOption(AZ_output, AZ_none);
solver.SetAztecOption(AZ_kspace, A->NumGlobalRows());
AztecOO_Operator AOpInv(&solver, A->NumGlobalRows());
Epetra_InvOperator AInvOp(&AOpInv);
EPETRA_CHK_ERR(EpetraExt::OperatorToMatlabFile("Ainv.dat", AInvOp));
comm.Barrier();
EPETRA_CHK_ERR(EpetraExt::MatlabFileToCrsMatrix("Ainv.dat", comm, AInv));
EPETRA_CHK_ERR(AInv->Apply(b,x));
EPETRA_CHK_ERR(x.Update(1.0, xx, -1.0));
double residual = 0.0;
EPETRA_CHK_ERR(x.Norm2(&residual));
if (verbose) cout << "Norm of difference between computed x and exact x = " << residual << endl;
int ierr = checkValues(residual,0.0,"Norm of difference between computed A1x1 and A1x1 from file", verbose);
delete A;
delete AInv;
#ifdef HAVE_MPI
MPI_Finalize() ;
#endif
return(ierr);
}

Definition at line 115 of file EpetraExt_MatlabEngine.h.

Constructor & Destructor Documentation

EpetraExt::EpetraExt_MatlabEngine::EpetraExt_MatlabEngine ( const Epetra_Comm Comm)

EpetraExt_MatlabEngine constructor which creates a MatlabEngine object with a connection to an instance of the application Matlab by starting a new Matlab process.

Parameters
Comm(In) An Epetra_Comm object.
Returns
A MatlabEngine object

Definition at line 65 of file EpetraExt_MatlabEngine.cpp.

EpetraExt::EpetraExt_MatlabEngine::~EpetraExt_MatlabEngine ( void  )

EpetraExt_MatlabEngine destructor which closes the connection to Matlab which causes the Matlab process to also exit.

Definition at line 73 of file EpetraExt_MatlabEngine.cpp.

Member Function Documentation

int EpetraExt::EpetraExt_MatlabEngine::EvalString ( char *  command,
char *  outputBuffer = NULL,
int  outputBufferSize = -1 
)

Sends a command to Matlab.

Any command that can normally be typed in at the Matlab command line can be passed in to EvalString(char* command). Commands such as 'help desk', 'edit', and 'plot(MATRIX)' will pop up an interactive window.

Parameters
command(In) the matlab command to run
outputBuffer(Out) (Optional) a user preallocated buffer for Matlab text output
outputBufferSize(In) (Optional) the size of the outputBuffer
Returns
Error Codes, see Detailed Description for more information

Definition at line 81 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutMultiVector ( const Epetra_MultiVector A,
const char *  variableName 
)

Puts a copy of the serial or distributed MultiVector into the Matlab workspace.

Parameters
A(In) the Epetra_MultiVector to put into Matlab
variableName(In) the variable name in the Matlab workspace of the Matlab double array (matrix) that will contain the values of the MultiVector
Returns
Error Codes, see Detailed Description for more information

Definition at line 98 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutRowMatrix ( const Epetra_RowMatrix A,
const char *  variableName,
bool  transA 
)

Puts a copy of the serial or distributed RowMatrix into the Matlab workspace.

Parameters
A(In) the Epetra_RowMatrix to put into Matlab
variableName(In) the variable name in the Matlab workspace of the Matlab sparse double array (matrix) that will contain the values of the RowMatrix
transA(In) if true then the transpose of A is put into Matlab NOTE: It is faster to put the transpose of A into Matlab since Matlab stores matrices in column-major form whereas Epetra stores them in row-major form.
Returns
Error Codes, see Detailed Description for more information

Definition at line 122 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutCrsGraph ( const Epetra_CrsGraph A,
const char *  variableName,
bool  transA 
)

not implemented yet

Definition at line 170 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutSerialDenseMatrix ( const Epetra_SerialDenseMatrix A,
const char *  variableName,
int  proc = 0 
)

Puts a copy of the SerialDenseMatrix into the Matlab workspace.

Parameters
A(In) the Epetra_SerialDenseMatrix to put into Matlab
variableName(In) the variable name in the Matlab workspace of the Matlab double array (matrix) that will contain the values of the SerialDenseMatrix
proc(In) for serial environment set to 0 for a parallel environment set to the process ID that owns the SerialDenseMatrix
Warning
The same parameters must be passed to each process.
Returns
Error Codes, see Detailed Description for more information

Definition at line 175 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutIntSerialDenseMatrix ( const Epetra_IntSerialDenseMatrix A,
const char *  variableName,
int  proc = 0 
)

Puts a copy of the IntSerialDenseMatrix into the Matlab workspace.

Parameters
A(In) the Epetra_IntSerialDenseMatrix to put into Matlab
variableName(In) the variable name in the Matlab workspace of the Matlab double array (matrix) that will contain the values of the IntSerialDenseMatrix
proc(In) for serial environment set to 0 for a parallel environment set to the process ID that owns the IntSerialDenseMatrix
Warning
The same parameters must be passed to each process.
Returns
Error Codes, see Detailed Description for more information

Definition at line 235 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutBlockMap ( const Epetra_BlockMap blockMap,
const char *  variableName,
bool  transA 
)

Puts a copy of the BlockMap or Map into the Matlab workspace.

Parameters
blockMap(In) the Epetra_BlockMap to put into Matlab
variableName(In) the variable name in the Matlab workspace of the Matlab sparse double array (matrix) that will contain the values of the BlockMap
transA(In) if true then the transpose of blockMap is put into Matlab NOTE: It is faster to put the transpose of blockMap into Matlab since Matlab stores matrices in column-major form whereas Epetra stores them in row-major form.
Returns
Error Codes, see Detailed Description for more information

Definition at line 339 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::PutIntoMatlab ( const char *  variableName,
mxArray *  matlabA 
)

Puts a mxArray into Matlab.

The Matlab provided C library provides mxArray which is used to construct and fill a Matlab object before sending it to Matlab to be put into the Matlab workspace. The mxArray is copied into the Matlab memory space and is not needed after it has been passed to Matlab. mxArrays should be destroyed using mxDestoryArray(mxArray) when they are no longer needed by the C/C++ program using them. Objects in Matlab must be destroyed using EvalString(char* command) and the appropriate Matlab command to destroy the object. EpetraExt::MatlabEngine uses PutIntoMatlab to pass all mxArrays it generates into Matlab. However, a user can create, fill, and put his/her own mxArrays into Matlab using this method. To create a mxArray mex.h must be included. For more information on how to use mxArrays see Matlab's documentation (type helpdesk at the Matlab command prompt) and see the External API Reference section.

Parameters
variableName(In) the name for the mxArray once it has been put into the Matlab workspace
matlabA(In) the mxArray to put into the Matlab workspace
Returns
Matlab error code from engPutVariable for Matlab versions >= 6.5 or from engPutArray for all other versions

Definition at line 379 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::GetMultiVector ( const char *  variableName,
Epetra_MultiVector A 
)

Puts a Matlab variable into a MultiVector.

The values from the Matlab variable are exported to the MultiVector using an export object. Therefore the MultiVector must be prepared by the user just like any MultiVector would be before calling an export.

Parameters
variableName(In) the name of the Matlab variable to be put into the given MultiVector
A(In) the MultiVector to put the Matlab variable values into
Returns
Error Codes, see Detailed Description for more information

Definition at line 394 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::GetSerialDenseMatrix ( const char *  variableName,
Epetra_SerialDenseMatrix A,
int  proc = 0 
)

Puts a Matlab variable into a SerialDenseMatrix on the specified PE.

The SerialDenseMatrix must be constructed by the user and have the proper amount of space to hold the values that will be copied from the given Matlab variable.

Parameters
variableName(In) the name of the Matlab variable to be put into the given SerialDenseMatrix
A(In) the SerialDenseMatrix to put the Matlab variable values into
proc(In) the PE that will own the SerialDenseMatrix
Returns
Error Codes, see Detailed Description for more information

Definition at line 435 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::GetIntSerialDenseMatrix ( const char *  variableName,
Epetra_IntSerialDenseMatrix A,
int  proc = 0 
)

Puts a Matlab variable into a IntSerialDenseMatrix on the specified PE.

The IntSerialDenseMatrix must be constructed by the user and have the proper amount of space to hold the values that will be copied from the given Matlab variable.

Parameters
variableName(In) the name of the Matlab variable to be put into the given IntSerialDenseMatrix
A(In) the IntSerialDenseMatrix to put the Matlab variable values into
proc(In) the PE that will own the IntSerialDenseMatrix
Returns
Error Codes, see Detailed Description for more information

Definition at line 460 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::GetCrsMatrix ( const char *  variableName,
Epetra_CrsMatrix A,
bool  getTrans 
)

Puts a Matlab variable into a CrsMatrix.

The values from the Matlab variable are exported to the CrsMatrix using an export object. Therefore the CrsMatrix must be prepared by the user just like any CrsMatrix would be before calling an export.

Warning
Getting a CrsMatrix from Matlab will cause a temporary Matlab variable to be created and then deleted named TRANS_variableName where variableName is the name of the Matlab variable.
Parameters
variableName(In) the name of the Matlab variable to be put into the given CrsMatrix
A(In) the CrsMatrix to put the Matlab variable values into
getTrans(In) if false then a temporary Matlab variable is created TRANS_variableName if true then the transpose of A copied from Matlab NOTE: It is faster to copy the transpose of A from Matlab since Matlab stores matrices in column-major form whereas Epetra stores them in row-major form.
Returns
Error Codes, see Detailed Description for more information

Definition at line 491 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::GetmxArrayDimensions ( mxArray *  matlabA,
bool &  isSparse,
int &  numRows,
int &  numCols,
int &  numNonZeros 
)

Get general information about the mxArray. For internal use but can be used by an advanced user.

Makes several Matlab function calls on the mxArray in order to determine the number of rows, columns, nonzeros, and whether or not the mxArray is sparse.

Parameters
matlabA(In) the mxArray to get information about
isSparse(Out) true if the mxArray is sparse
numRows(Out) the number of rows in the mxArray
numCols(Out) the number of columns in the mxArray
numNonZeros(Out) the number of nonzeros in the mxArray
Returns
Error Codes, see Detailed Description for more information

Definition at line 576 of file EpetraExt_MatlabEngine.cpp.

int EpetraExt::EpetraExt_MatlabEngine::GetmxArray ( const char *  variableName,
mxArray **  matlabA 
)

Get a mxArray from Matlab. For internal use but can be used by an advanced user.

Calls the Matlab provided engGetVariable function which copies the specified variable from the Matlab workspace into a mxArray. Hence any changes to mxArray will not show up in Matlab, and any changes to the Matlab variable in Matlab will not show up in the mxArray. When finished with the mxArray mxDestroyArray should be called on the mxArray. Matlab appears to perform a copy to a new mxArray each time engGetVariable is called.

Parameters
variableName(In) the mxArray to get from the Matlab workspace
matlabA(Out) a pointer to a mxArray* to put the mxArray that engGetVariable returns
Returns
Error Codes, see Detailed Description for more information

Definition at line 599 of file EpetraExt_MatlabEngine.cpp.


The documentation for this class was generated from the following files: