TrilinosCouplings  Development
 All Classes Namespaces Files Functions Pages
Classes | Macros | Typedefs | Functions
example_Maxwell.cpp File Reference

Example solution of the eddy current Maxwell's equations using curl-conforming (edge) elements. More...

#include "TrilinosCouplings_config.h"
#include "TrilinosCouplings_Pamgen_Utils.hpp"
#include "TrilinosCouplings_Statistics.hpp"
#include "Intrepid_FunctionSpaceTools.hpp"
#include "Intrepid_FieldContainer.hpp"
#include "Intrepid_CellTools.hpp"
#include "Intrepid_ArrayTools.hpp"
#include "Intrepid_HCURL_HEX_I1_FEM.hpp"
#include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
#include "Intrepid_RealSpaceTools.hpp"
#include "Intrepid_DefaultCubatureFactory.hpp"
#include "Intrepid_Utils.hpp"
#include "Kokkos_Core.hpp"
#include "Epetra_Time.h"
#include "Epetra_Map.h"
#include "Epetra_SerialComm.h"
#include "Epetra_Import.h"
#include "Epetra_Export.h"
#include "Epetra_CrsMatrix.h"
#include "Epetra_FECrsMatrix.h"
#include "Epetra_FEVector.h"
#include "Epetra_LinearProblem.h"
#include "Teuchos_oblackholestream.hpp"
#include "Teuchos_RCP.hpp"
#include "Teuchos_BLAS.hpp"
#include "Teuchos_GlobalMPISession.hpp"
#include "Teuchos_ParameterList.hpp"
#include "Teuchos_XMLParameterListHelpers.hpp"
#include "Teuchos_Comm.hpp"
#include "Teuchos_DefaultComm.hpp"
#include "Shards_CellTopology.hpp"
#include "EpetraExt_RowMatrixOut.h"
#include "EpetraExt_MultiVectorOut.h"
#include "EpetraExt_VectorOut.h"
#include "EpetraExt_MatrixMatrix.h"
#include "AztecOO.h"
#include "ml_MultiLevelPreconditioner.h"
#include "ml_RefMaxwell_11_Operator.h"
#include "ml_RefMaxwell.h"
#include "ml_EdgeMatrixFreePreconditioner.h"
#include "ml_epetra_utils.h"
#include "create_inline_mesh.h"
#include "pamgen_im_exodusII_l.h"
#include "pamgen_im_ne_nemesisI_l.h"
#include "pamgen_extras.h"
#include "TrilinosCouplings_IntrepidPoissonExampleHelpers.hpp"
Include dependency graph for example_Maxwell.cpp:

Classes

struct  fecomp
 

Macros

#define ABS(x)   ((x)>0?(x):-(x))
 
#define SQR(x)   ((x)*(x))
 

Typedefs

typedef
Intrepid::FunctionSpaceTools 
IntrepidFSTools
 
typedef
Intrepid::RealSpaceTools
< double > 
IntrepidRSTools
 
typedef Intrepid::CellTools
< double > 
IntrepidCTools
 

Functions

template<class Container >
double distance (Container &nodeCoord, int i1, int i2)
 
void TestMultiLevelPreconditioner_Maxwell (char ProblemType[], Teuchos::ParameterList &MLList, Epetra_CrsMatrix &CurlCurl, Epetra_CrsMatrix &D0clean, Epetra_CrsMatrix &M0inv, Epetra_CrsMatrix &Ms, Epetra_CrsMatrix &M1, Epetra_MultiVector &xh, Epetra_MultiVector &b, double &TotalErrorResidual, double &TotalErrorExactSol)
 ML Preconditioner. More...
 
int evalu (double &uExact0, double &uExact1, double &uExact2, double &x, double &y, double &z)
 ML Preconditioner. More...
 
int evalCurlu (double &curlu0, double &curlu1, double &curlu2, double &x, double &y, double &z, double &mu)
 Curl of exact solution. More...
 
int evalCurlCurlu (double &curlcurlu0, double &curlcurlu1, double &curlcurlu2, double &x, double &y, double &z, double &mu)
 CurlCurl of exact solution. More...
 
int body (int argc, char *argv[])
 
int main (int argc, char *argv[])
 
int Multiply_Ones (const Epetra_CrsMatrix &A, const Epetra_Vector &x, Epetra_Vector &y)
 Multiplies Ax = y, where all non-zero entries of A are replaced with the value 1.0. More...
 
void solution_test (string msg, const Epetra_Operator &A, const Epetra_MultiVector &lhs, const Epetra_MultiVector &rhs, const Epetra_MultiVector &xexact, Epetra_Time &Time, double &TotalErrorExactSol, double &TotalErrorResidual)
 Compute ML solution residual. More...
 
void TestMueLuMultiLevelPreconditioner_Maxwell (char ProblemType[], Teuchos::ParameterList &MLList, Epetra_CrsMatrix &CurlCurl, Epetra_CrsMatrix &D0clean, Epetra_CrsMatrix &M0inv, Epetra_CrsMatrix &Ms, Epetra_CrsMatrix &M1, Epetra_MultiVector &coords, Epetra_MultiVector &xh, Epetra_MultiVector &b, double &TotalErrorResidual, double &TotalErrorExactSol)
 

Detailed Description

Example solution of the eddy current Maxwell's equations using curl-conforming (edge) elements.

This example uses the following Trilinos packages:

Maxwell System:

curl x mu^{-1} curl E + sigma E = f

Corresponding discrete linear system for edge element coeficients (x):

(Kc + Mc)x = b

Kc    - Hcurl stiffness matrix
Mc    - Hcurl mass matrix
b     - right hand side vector
Author
Created by P. Bochev, D. Ridzal, K. Peterson, D. Hensinger, C. Siefert.
Remarks
Usage
./TrilinosCouplings_examples_scaling_Example_Maxwell.exe  inputfile.xml


inputfile.xml (optional)  -  xml input file containing Pamgen mesh description
and material parameters for each Pamgen block,
if not present code attempts to read Maxwell.xml.

Input files available in Trilinos for use with the Maxwell driver:

Function Documentation

int evalCurlCurlu ( double &  curlcurlu0,
double &  curlcurlu1,
double &  curlcurlu2,
double &  x,
double &  y,
double &  z,
double &  mu 
)

CurlCurl of exact solution.

Parameters
curlcurlu0[out] first component of curl-curl of exact solution
curlcurlu1[out] second component of curl-curl of exact solution
curlcurlu2[out] third component of curl-curl of exact solution
x[in] x coordinate
y[in] y coordinate
z[in] z coordinate
mu[in] material parameter
int evalCurlu ( double &  curlu0,
double &  curlu1,
double &  curlu2,
double &  x,
double &  y,
double &  z,
double &  mu 
)

Curl of exact solution.

Parameters
curlu0[out] first component of curl of exact solution
curlu1[out] second component of curl of exact solution
curlu2[out] third component of curl of exact solution
x[in] x coordinate
y[in] y coordinate
z[in] z coordinate
mu[in] material parameter
int evalu ( double &  uExact0,
double &  uExact1,
double &  uExact2,
double &  x,
double &  y,
double &  z 
)

ML Preconditioner.

Parameters
ProblemType[in] problem type
MLList[in] ML parameter list
CurlCurl[in] H(curl) stiffness matrix
D0clean[in] Edge to node stiffness matrix
M0inv[in] H(grad) mass matrix inverse
Ms[in] H(curl) mass matrix w/ sigma
M1[in] H(curl) mass matrix w/o sigma
xh[out] solution vector
b[in] right-hand-side vector
TotalErrorResidual[out] error residual
TotalErrorExactSol[out] error in xhExact solution evaluation.
uExact0[out] first component of exact solution at (x,y,z)
uExact1[out] second component of exact solution at (x,y,z)
uExact2[out] third component of exact solution at (x,y,z)
x[in] x coordinate
y[in] y coordinate
z[in] z coordinate
int Multiply_Ones ( const Epetra_CrsMatrix &  A,
const Epetra_Vector &  x,
Epetra_Vector &  y 
)

Multiplies Ax = y, where all non-zero entries of A are replaced with the value 1.0.

Parameters
A[in] matrix
x[in] vector
y[in] vector
void solution_test ( string  msg,
const Epetra_Operator &  A,
const Epetra_MultiVector &  lhs,
const Epetra_MultiVector &  rhs,
const Epetra_MultiVector &  xexact,
Epetra_Time &  Time,
double &  TotalErrorExactSol,
double &  TotalErrorResidual 
)

Compute ML solution residual.

Parameters
A[in] discrete operator
lhs[in] solution vector
rhs[in] right hand side vector
Time[in] elapsed time for output
TotalErrorResidual[out] error residual
TotalErrorExactSol[out] error in xh (not an appropriate measure for H(curl) basis functions)
void TestMultiLevelPreconditioner_Maxwell ( char  ProblemType[],
Teuchos::ParameterList &  MLList,
Epetra_CrsMatrix &  CurlCurl,
Epetra_CrsMatrix &  D0clean,
Epetra_CrsMatrix &  M0inv,
Epetra_CrsMatrix &  Ms,
Epetra_CrsMatrix &  M1,
Epetra_MultiVector &  xh,
Epetra_MultiVector &  b,
double &  TotalErrorResidual,
double &  TotalErrorExactSol 
)

ML Preconditioner.

Parameters
ProblemType[in] problem type
MLList[in] ML parameter list
CurlCurl[in] H(curl) stiffness matrix
D0clean[in] Edge to node stiffness matrix
M0inv[in] H(grad) mass matrix inverse
Ms[in] H(curl) mass matrix w/ sigma
M1[in] H(curl) mass matrix w/o sigma
xh[out] solution vector
b[in] right-hand-side vector
TotalErrorResidual[out] error residual
TotalErrorExactSol[out] error in xh

References solution_test().