TrilinosCouplings Development
Loading...
Searching...
No Matches
example_Maxwell_Tpetra.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 "TrilinosCouplings_IntrepidPoissonExample_SolveWithBelos.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 "Tpetra_Map.hpp"
#include "Tpetra_Import.hpp"
#include "Tpetra_Export.hpp"
#include "Tpetra_CrsMatrix.hpp"
#include "Tpetra_FECrsMatrix.hpp"
#include "Tpetra_FECrsGraph.hpp"
#include "Tpetra_FEMultiVector.hpp"
#include "Tpetra_Assembly_Helpers.hpp"
#include "TpetraExt_MatrixMatrix.hpp"
#include "MatrixMarket_Tpetra.hpp"
#include "Tpetra_applyDirichletBoundaryCondition.hpp"
#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_OrdinalTraits.hpp"
#include "Teuchos_StackedTimer.hpp"
#include "Shards_CellTopology.hpp"
#include <Xpetra_MultiVector.hpp>
#include <Xpetra_MultiVectorFactory.hpp>
#include <Xpetra_CrsMatrix.hpp>
#include <Xpetra_CrsMatrixWrap.hpp>
#include <Xpetra_Matrix.hpp>
#include <BelosConfigDefs.hpp>
#include <BelosLinearProblem.hpp>
#include <BelosSolverFactory.hpp>
#include <BelosXpetraAdapter.hpp>
#include <BelosMueLuAdapter.hpp>
#include <MueLu_RefMaxwell.hpp>
#include <MueLu_Exceptions.hpp>
#include "create_inline_mesh.h"
#include "pamgen_im_exodusII_l.h"
#include "pamgen_im_ne_nemesisI_l.h"
#include "pamgen_extras.h"
Include dependency graph for example_Maxwell_Tpetra.cpp:

Classes

struct  fecomp

Macros

#define ABS(x)
#define SQR(x)
#define TC_sumAll(rcpComm, in, out)
#define TC_minAll(rcpComm, in, out)
#define TC_maxAll(rcpComm, in, out)

Typedefs

typedef double SC
typedef int LO
typedef Tpetra::Map ::global_ordinal_type GO
typedef Tpetra::Map ::node_type Node
typedef Node NO
typedef Tpetra::CrsMatrix< SC, LO, GO, Node > Tpetra_CrsMatrix
typedef Tpetra::FECrsMatrix< SC, LO, GO, Node > Tpetra_FECrsMatrix
typedef Tpetra::FECrsGraph< LO, GO, Node > Tpetra_FECrsGraph
typedef Tpetra::MultiVector< SC, LO, GO, Node > Tpetra_MultiVector
typedef Tpetra::FEMultiVector< SC, LO, GO, Node > Tpetra_FEMultiVector
typedef Tpetra::Vector< SC, LO, GO, Node > Tpetra_Vector
typedef Tpetra::Map< LO, GO, Node > Tpetra_Map
typedef Tpetra::Import< LO, GO, Node > Tpetra_Import
typedef Tpetra::Export< LO, GO, Node > Tpetra_Export
typedef Tpetra::Operator< SC, LO, GO, Node > Tpetra_Operator
typedef Xpetra::Operator< SC, LO, GO, Node > Xpetra_Operator
typedef Xpetra::Matrix< SC, LO, GO, NO > Matrix
typedef Xpetra::MultiVector< SC, LO, GO, NO > Xpetra_MultiVector
typedef Intrepid::RealSpaceTools< double > IntrepidRSTools
typedef Intrepid::CellTools< double > IntrepidCTools

Functions

template<class Container>
double distance (Container &nodeCoord, int i1, int i2)
RCP< Matrix > toXpetra (RCP< Tpetra_CrsMatrix > &mat)
RCP< Xpetra_MultiVector > toXpetra (RCP< Tpetra_MultiVector > &vec)
RCP< Xpetra_Operator > BuildPreconditioner_MueLu (char ProblemType[], Teuchos::ParameterList &MLList, RCP< Tpetra_CrsMatrix > &CurlCurl, RCP< Tpetra_CrsMatrix > &D0clean, RCP< Tpetra_CrsMatrix > &M0inv, RCP< Tpetra_CrsMatrix > &Ms, RCP< Tpetra_CrsMatrix > &M1, RCP< Tpetra_MultiVector > &coords)
 MueLu Preconditioner.
int evalu (double &uExact0, double &uExact1, double &uExact2, double &x, double &y, double &z)
 Exact solution evaluation.
int evalCurlu (double &curlu0, double &curlu1, double &curlu2, double &x, double &y, double &z, double &mu)
 Curl of exact solution.
int evalCurlCurlu (double &curlcurlu0, double &curlcurlu1, double &curlcurlu2, double &x, double &y, double &z, double &mu)
 CurlCurl of exact solution.
int main (int argc, char *argv[])
void solution_test (string msg, const Tpetra_Operator &A, const Tpetra_MultiVector &lhs, const Tpetra_MultiVector &rhs, const Tpetra_MultiVector &xexact, double &TotalErrorExactSol, double &TotalErrorResidual)
 Compute ML solution residual.

Detailed Description

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

This example uses the following Trilinos packages:

  • Pamgen to generate a Hexahedral mesh.
  • Intrepid to build the discretization matrices and right-hand side.
  • Epetra to handle the global matrix and vector.
  • ML to solve the linear system.
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:

  • Maxwell.xml - basic input file with box mesh and one mesh block
  • Ninja.xml - input file with distorted mesh (shaped like Ninja star) and one mesh block
  • Maxwell_in_block.xml - input file with box mesh with a center block with different material values.

Macro Definition Documentation

◆ ABS

#define ABS ( x)
Value:
((x)>0?(x):-(x))

◆ SQR

#define SQR ( x)
Value:
((x)*(x))

◆ TC_maxAll

#define TC_maxAll ( rcpComm,
in,
out )
Value:
Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MAX, in, Teuchos::outArg(out))

◆ TC_minAll

#define TC_minAll ( rcpComm,
in,
out )
Value:
Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_MIN, in, Teuchos::outArg(out))

◆ TC_sumAll

#define TC_sumAll ( rcpComm,
in,
out )
Value:
Teuchos::reduceAll(*rcpComm, Teuchos::REDUCE_SUM, in, Teuchos::outArg(out))

Function Documentation

◆ BuildPreconditioner_MueLu()

RCP< Xpetra_Operator > BuildPreconditioner_MueLu ( char ProblemType[],
Teuchos::ParameterList & MLList,
RCP< Tpetra_CrsMatrix > & CurlCurl,
RCP< Tpetra_CrsMatrix > & D0clean,
RCP< Tpetra_CrsMatrix > & M0inv,
RCP< Tpetra_CrsMatrix > & Ms,
RCP< Tpetra_CrsMatrix > & M1,
RCP< Tpetra_MultiVector > & coords )

MueLu Preconditioner.

Parameters
ProblemType[in] problem type
MLList[in] 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 sigm

◆ evalCurlCurlu()

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

◆ evalCurlu()

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

◆ evalu()

int evalu ( double & uExact0,
double & uExact1,
double & uExact2,
double & x,
double & y,
double & z )

Exact solution evaluation.

Parameters
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

◆ solution_test()

void solution_test ( string msg,
const Tpetra_Operator & A,
const Tpetra_MultiVector & lhs,
const Tpetra_MultiVector & rhs,
const Tpetra_MultiVector & xexact,
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)