Example solution of a div-curl system on a hexahedral mesh using div-conforming (face) elements. More...
#include "TrilinosCouplings_config.h"#include "TrilinosCouplings_Pamgen_Utils.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_HDIV_HEX_I1_FEM.hpp"#include "Intrepid_RealSpaceTools.hpp"#include "Intrepid_DefaultCubatureFactory.hpp"#include "Intrepid_Utils.hpp"#include "Epetra_Time.h"#include "Epetra_Map.h"#include "Epetra_SerialComm.h"#include "Epetra_FECrsMatrix.h"#include "Epetra_FEVector.h"#include "Epetra_Vector.h"#include "Epetra_LinearProblem.h"#include "Epetra_Import.h"#include "Epetra_Export.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 "Shards_CellTopology.hpp"#include "EpetraExt_RowMatrixOut.h"#include "EpetraExt_MultiVectorOut.h"#include "EpetraExt_VectorOut.h"#include "EpetraExt_MatrixMatrix.h"#include "create_inline_mesh.h"#include "pamgen_im_exodusII_l.h"#include "pamgen_im_ne_nemesisI_l.h"#include "pamgen_extras.h"#include "AztecOO.h"#include "ml_epetra_utils.h"#include "ml_RefMaxwell_11_Operator.h"#include "ml_FaceMatrixFreePreconditioner.h"#include "ml_RefMaxwell.h"
Classes | |
| struct | fecomp |
Macros | |
| #define | ABS(x) ((x)>0?(x):-(x)) |
Typedefs | |
|
typedef Intrepid::FunctionSpaceTools | IntrepidFSTools |
|
typedef Intrepid::RealSpaceTools < double > | IntrepidRSTools |
|
typedef Intrepid::CellTools < double > | IntrepidCTools |
Functions | |
| int | Multiply_Abs (const Epetra_CrsMatrix &A, const Epetra_Vector &x, Epetra_Vector &y) |
| Multiplies abs(A)x = y, where all non-zero entries of A are replaced with their absolute values value. More... | |
| void | TestMultiLevelPreconditioner_DivLSFEM (char ProblemType[], Teuchos::ParameterList &MLList, Epetra_CrsMatrix &GradDiv, Epetra_CrsMatrix &D0clean, Epetra_CrsMatrix &D1clean, Epetra_CrsMatrix &FaceNode, Epetra_CrsMatrix &M1, Epetra_CrsMatrix &M1inv, Epetra_CrsMatrix &M2, 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) |
| Exact solution evaluation. More... | |
| double | evalDivu (double &x, double &y, double &z) |
| Divergence of exact solution. 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) |
| Curl of curl of exact solution. More... | |
| 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... | |
Example solution of a div-curl system on a hexahedral mesh using div-conforming (face) elements.
This example uses Pamgen to generate a hexahedral mesh, Intrepid to
build mass and stiffness matrices, and ML to solve.
System
div v + \phi = f in Omega
v+ A grad \phi = 0 in Omega
Dirichlet BC: \phi given on Gamma
Omega is the box [0,1]^3
where f is derived from a prescribed solution
\phi(x,y,z)=-sin^2(\pi x)*sin^2(\pi y)*sin^2(\pi z)
or other, discontinuous solution, aka "patch test" ./TrilinosCouplings_examples_scaling_Example_DivLSFEM.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 DivLSFEMin.xml.This example uses the following Trilinos packages:
For more details on the formulation see
Div-Curl System:
curl u = g in Omega
div u = h in Omega
u.n = 0 on Gamma
Corresponding discrete linear system for face element coeficients (x):
(Kd + Md*Dc*McInv*Dc'*Md)x = b
Kd - Hdiv stiffness matrix
Md - Hdiv mass matrix
Dc - Edge to Face incidence matrix
McInv - Hgrad mass matrix inverse
b - right hand side vector
./TrilinosCouplings_examples_scaling_Example_DivLSFEM.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 DivLSFEMin.xml. | int evalCurlCurlu | ( | double & | curlCurlu0, |
| double & | curlCurlu1, | ||
| double & | curlCurlu2, | ||
| double & | x, | ||
| double & | y, | ||
| double & | z, | ||
| double & | mu | ||
| ) |
Curl of curl of exact solution.
| 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.
| 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 |
| double evalDivu | ( | double & | x, |
| double & | y, | ||
| double & | z | ||
| ) |
Divergence of exact solution.
| x | [in] x coordinate |
| y | [in] y coordinate |
| z | [in] z coordinate |
| int evalu | ( | double & | uExact0, |
| double & | uExact1, | ||
| double & | uExact2, | ||
| double & | x, | ||
| double & | y, | ||
| double & | z | ||
| ) |
Exact 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_Abs | ( | const Epetra_CrsMatrix & | A, |
| const Epetra_Vector & | x, | ||
| Epetra_Vector & | y | ||
| ) |
Multiplies abs(A)x = y, where all non-zero entries of A are replaced with their absolute values value.
| A | [in] matrix |
| x | [in] vector |
| y | [out] vector |
| 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.
| 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.
| 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(div) basis functions) |
| void TestMultiLevelPreconditioner_DivLSFEM | ( | char | ProblemType[], |
| Teuchos::ParameterList & | MLList, | ||
| Epetra_CrsMatrix & | GradDiv, | ||
| Epetra_CrsMatrix & | D0clean, | ||
| Epetra_CrsMatrix & | D1clean, | ||
| Epetra_CrsMatrix & | FaceNode, | ||
| Epetra_CrsMatrix & | M1, | ||
| Epetra_CrsMatrix & | M1inv, | ||
| Epetra_CrsMatrix & | M2, | ||
| Epetra_MultiVector & | xh, | ||
| Epetra_MultiVector & | b, | ||
| double & | TotalErrorResidual, | ||
| double & | TotalErrorExactSol | ||
| ) |
ML Preconditioner.
| ProblemType | [in] problem type |
| MLList | [in] ML parameter list |
| GradDiv | [in] H(curl) stiffness matrix |
| D0clean | [in] Edge to node incidence matrix |
| D1clean | [in] Face to edge incidence matrix |
| FaceNode | [in] Face to node incidence matrix |
| M1inv | [in] H(curl) mass matrix inverse |
| M2 | [in] H(div) mass matrix |
| xh | [out] solution vector |
| b | [in] right-hand-side vector |
| TotalErrorResidual | [out] error residual |
| TotalErrorExactSol | [out] error in xh |
References Multiply_Ones(), and solution_test().
1.8.5