TrilinosCouplings  Development
 All Classes Namespaces Files Functions Pages
Typedefs | Functions | Variables
example_StabilizedADR.cpp File Reference

Example solution of a steady-state advection-diffusion-reaction equation with Dirichlet boundary conditon on a hexahedral mesh using nodal (Hgrad) elements and stabilization. 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_HGRAD_HEX_C1_FEM.hpp"
#include "Intrepid_RealSpaceTools.hpp"
#include "Intrepid_DefaultCubatureFactory.hpp"
#include "Intrepid_Utils.hpp"
#include "Sacado_No_Kokkos.hpp"
#include "Epetra_Time.h"
#include "Epetra_Map.h"
#include "Epetra_SerialComm.h"
#include "Epetra_FECrsMatrix.h"
#include "Epetra_FEVector.h"
#include "Epetra_Import.h"
#include "Teuchos_oblackholestream.hpp"
#include "Teuchos_RCP.hpp"
#include "Teuchos_BLAS.hpp"
#include "Teuchos_GlobalMPISession.hpp"
#include "Teuchos_XMLParameterListHelpers.hpp"
#include "Teuchos_Assert.hpp"
#include "Shards_CellTopology.hpp"
#include "EpetraExt_RowMatrixOut.h"
#include "EpetraExt_MultiVectorOut.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_MultiLevelPreconditioner.h"
#include "ml_epetra_utils.h"
#include "map"
Include dependency graph for example_StabilizedADR.cpp:

Typedefs

typedef Sacado::Fad::SFad
< double, 3 > 
Fad3
 
typedef shards::CellTopology ShardsCellTopology
 
typedef
Intrepid::FunctionSpaceTools 
IntrepidFSTools
 
typedef Intrepid::CellTools
< double > 
IntrepidCTools
 
typedef long long int128
 

Functions

template<typename Scalar >
const Scalar exactSolution (const Scalar &x, const Scalar &y, const Scalar &z)
 User-defined exact solution. More...
 
template<typename Scalar >
void diffusionTensor (Scalar diffusion[][3], const Scalar &x, const Scalar &y, const Scalar &z)
 User-defined diffusion tensor. More...
 
template<typename Scalar >
void advectiveVector (Scalar advection[3], const Scalar &x, const Scalar &y, const Scalar &z)
 User-defined advective vector. More...
 
template<typename Scalar >
const Scalar reactionTerm (const Scalar &x, const Scalar &y, const Scalar &z)
 User-defined reaction coefficient. More...
 
template<typename Scalar >
void exactSolutionGrad (Scalar gradExact[3], const Scalar &x, const Scalar &y, const Scalar &z)
 Computes gradient of the exact solution. Requires user-defined exact solution. More...
 
template<typename Scalar >
const Scalar sourceTerm (Scalar &x, Scalar &y, Scalar &z)
 Computes source term: f = -div(K.grad u - b.u) + c.u. Requires user-defined exact solution, diffusion tensor, advective vector and reaction coefficient. More...
 
template<class ArrayOut , class ArrayIn >
void evaluateDiffusionTensor (ArrayOut &worksetDiffusionValues, const ArrayIn &evaluationPoints)
 Computation of the diffusion tensor at array of points in physical space. More...
 
template<class ArrayOut , class ArrayIn >
void evaluateAdvectiveVector (ArrayOut &worksetAdvectionValues, const ArrayIn &evaluationPoints)
 Computation of the advective vector at array of points in physical space.. More...
 
template<class ArrayOut , class ArrayIn >
void evaluateReactionCoefficient (ArrayOut &worksetReactionValues, const ArrayIn &evaluationPoints)
 Computation of the reaction coefficient at array of points in physical space.. More...
 
template<class ArrayOut , class ArrayIn >
void evaluateSourceTerm (ArrayOut &sourceTermValues, const ArrayIn &evaluationPoints)
 Computation of the source term at array of points in physical space. More...
 
template<class ArrayOut , class ArrayIn >
void evaluateExactSolution (ArrayOut &exactSolutionValues, const ArrayIn &evaluationPoints)
 Computation of the exact solution at array of points in physical space. More...
 
template<class ArrayOut , class ArrayIn >
void evaluateExactSolutionGrad (ArrayOut &exactSolutionGradValues, const ArrayIn &evaluationPoints)
 Computation of the gradient of the exact solution at array of points in physical space. More...
 
int TestMultiLevelPreconditioner (char ProblemType[], Teuchos::ParameterList &MLList, Epetra_CrsMatrix &A, const Epetra_MultiVector &xexact, Epetra_MultiVector &b, Epetra_MultiVector &uh, double &TotalErrorResidual, double &TotalErrorExactSol)
 
template<class Scalar >
void getPamgenMesh (FieldContainer< Scalar > &localNodeCoordsFC, FieldContainer< long long > &localCellToNodeFC, FieldContainer< int > &nodeOnBoundaryFC, FieldContainer< bool > &nodeIsOwnedFC, FieldContainer< long long > &globalNodeIdsFC, const std::string &meshInput, const int &procRank, const int &numProcs, const Epetra_Comm &Comm, Epetra_Time &Time, const string &message, const int verbose=0)
 Pamgen wrapper. More...
 
void getInputArguments (Teuchos::ParameterList &inputMeshList, std::string &meshInput, Teuchos::ParameterList &inputSolverList, int &argc, char *argv[], const Epetra_Comm &Comm, const std::string &inputMeshFile="ADR.xml", const std::string &intputSolverFile="ML_nonsym.xml")
 
int main (int argc, char *argv[])
 

Variables

double g_advection [3] ={0.,0.,0.}
 
double g_reaction
 

Detailed Description

Example solution of a steady-state advection-diffusion-reaction equation with Dirichlet boundary conditon on a hexahedral mesh using nodal (Hgrad) elements and stabilization.

    This example requires the following Trilinos packages:
 Steady-state advection-diffusion boundary value problem in conservative form:

        -div (A.grad(u) - b.u) + c.u = f  in Omega
                                   u = 0  on Gamma
 where
        A   is a symmetric, positive definite diffusivity tensor
        b   is advective velocity vector which is not required to be divergence free
        c   is a reaction coefficient
        f   is a given source term

 The user must provide definitions for these quantities and definition for the exact solution


 Corresponding discrete linear system for nodal coefficients(x):

             Kx = d

        K - HGrad stiffness matrix
        d - right hand side vector
Author
Created by P. Bochev, D. Ridzal, K. Peterson, D. Hensinger, C. Siefert.
Remarks
Usage:
./TrilinosCouplings_examples_scaling_example_StabilizedADR.exe
Example requires Pamgen formatted mesh input file named ADR.xml.

Function Documentation

template<typename Scalar >
void advectiveVector ( Scalar  advection[3],
const Scalar &  x,
const Scalar &  y,
const Scalar &  z 
)

User-defined advective vector.

Parameters
advection[out] 3-dimensional advective vector evaluated at (x,y,z)
x[in] x-coordinate of the evaluation point
y[in] y-coordinate of the evaluation point
z[in] z-coordinate of the evaluation point
Remarks
Advective vector is not required to be divergence-free
template<typename Scalar >
void diffusionTensor ( Scalar  diffusion[][3],
const Scalar &  x,
const Scalar &  y,
const Scalar &  z 
)

User-defined diffusion tensor.

Parameters
diffusion[out] 3 x 3 diffusivity tensor evaluated at (x,y,z)
x[in] x-coordinate of the evaluation point
y[in] y-coordinate of the evaluation point
z[in] z-coordinate of the evaluation point
Warning
Symmetric and positive definite tensor is required for every (x,y,z).
template<class ArrayOut , class ArrayIn >
void evaluateAdvectiveVector ( ArrayOut &  worksetAdvectionValues,
const ArrayIn &  evaluationPoints 
)

Computation of the advective vector at array of points in physical space..

Parameters
worksetAdvectionValues[out] Rank-3 (C,P,D) array with the values of the advective vector
evaluationPoints[in] Rank-3 (C,P,D) array with the evaluation points in physical frame
template<class ArrayOut , class ArrayIn >
void evaluateDiffusionTensor ( ArrayOut &  worksetDiffusionValues,
const ArrayIn &  evaluationPoints 
)

Computation of the diffusion tensor at array of points in physical space.

Parameters
worksetDiffusionValues[out] Rank-2, 3 or 4 array with dimensions (C,P), (C,P,D) or (C,P,D,D) with the values of the diffusion tensor
evaluationPoints[in] Rank-3 (C,P,D) array with the evaluation points in physical frame
template<class ArrayOut , class ArrayIn >
void evaluateExactSolution ( ArrayOut &  exactSolutionValues,
const ArrayIn &  evaluationPoints 
)

Computation of the exact solution at array of points in physical space.

Parameters
exactSolutionValues[out] Rank-2 (C,P) array with the values of the exact solution
evaluationPoints[in] Rank-3 (C,P,D) array with the evaluation points in physical frame
template<class ArrayOut , class ArrayIn >
void evaluateExactSolutionGrad ( ArrayOut &  exactSolutionGradValues,
const ArrayIn &  evaluationPoints 
)

Computation of the gradient of the exact solution at array of points in physical space.

Parameters
exactSolutionGradValues[out] Rank-3 (C,P,D) array with the values of the gradient of the exact solution
evaluationPoints[in] Rank-3 (C,P,D) array with the evaluation points in physical frame
template<class ArrayOut , class ArrayIn >
void evaluateReactionCoefficient ( ArrayOut &  worksetReactionValues,
const ArrayIn &  evaluationPoints 
)

Computation of the reaction coefficient at array of points in physical space..

Parameters
worksetReactionValues[out] Rank-2 (C,P) array with the values of the reaction coefficient
evaluationPoints[in] Rank-3 (C,P,D) array with the evaluation points in physical frame
template<class ArrayOut , class ArrayIn >
void evaluateSourceTerm ( ArrayOut &  sourceTermValues,
const ArrayIn &  evaluationPoints 
)

Computation of the source term at array of points in physical space.

Parameters
sourceTermValues[out] Rank-2 (C,P) array with the values of the source term
evaluationPoints[in] Rank-3 (C,P,D) array with the evaluation points in physical frame
template<typename Scalar >
const Scalar exactSolution ( const Scalar &  x,
const Scalar &  y,
const Scalar &  z 
)

User-defined exact solution.

Parameters
x[in] x-coordinate of the evaluation point
y[in] y-coordinate of the evaluation point
z[in] z-coordinate of the evaluation point
Returns
Value of the exact solution at (x,y,z)
template<typename Scalar >
void exactSolutionGrad ( Scalar  gradExact[3],
const Scalar &  x,
const Scalar &  y,
const Scalar &  z 
)

Computes gradient of the exact solution. Requires user-defined exact solution.

Parameters
gradExact[out] gradient of the exact solution evaluated at (x,y,z)
x[in] x-coordinate of the evaluation point
y[in] y-coordinate of the evaluation point
z[in] z-coordinate of the evaluation point
template<class Scalar >
void getPamgenMesh ( FieldContainer< Scalar > &  localNodeCoordsFC,
FieldContainer< long long > &  localCellToNodeFC,
FieldContainer< int > &  nodeOnBoundaryFC,
FieldContainer< bool > &  nodeIsOwnedFC,
FieldContainer< long long > &  globalNodeIdsFC,
const std::string &  meshInput,
const int &  procRank,
const int &  numProcs,
const Epetra_Comm &  Comm,
Epetra_Time &  Time,
const string &  message,
const int  verbose = 0 
)

Pamgen wrapper.

Parameters
localNodeCoordsFC[out] Rank-2 (LN,D) array with the local nodes
localCellToNodeFC[out] Rank-2 (C,N) array with the local nodes connectivity
nodeOnBoundaryFC[out] Rank-1 (LN) array tells which local nodes are on the boundary
nodeIsOwnedFC[out] Rank-1 (LN) array tells which local nodes are owned
globalNodeIdsFC[out] Rank-1 (LN) array with the global node Ids of the local nodes
meshInput[in] The list of parameters for Create_Pamgen_Mesh, extracted from the mesh file
procRank[in] Rank of the processor
numProcs[in] Number of processors
Comm[in] Communicator object
Time[in] Time object used for timings of various mesh tasks in the wrapper
message[in] String that is printed as part of the various diagnostic outputs
verbose[in] Forces diagnostic output if set to 1, default is 0
template<typename Scalar >
const Scalar reactionTerm ( const Scalar &  x,
const Scalar &  y,
const Scalar &  z 
)

User-defined reaction coefficient.

Parameters
x[in] x-coordinate of the evaluation point
y[in] y-coordinate of the evaluation point
z[in] z-coordinate of the evaluation point
Returns
Value of the reaction coefficient function at (x,y,z)
template<typename Scalar >
const Scalar sourceTerm ( Scalar &  x,
Scalar &  y,
Scalar &  z 
)

Computes source term: f = -div(K.grad u - b.u) + c.u. Requires user-defined exact solution, diffusion tensor, advective vector and reaction coefficient.

Parameters
x[in] x-coordinate of the evaluation point
y[in] y-coordinate of the evaluation point
z[in] z-coordinate of the evaluation point
Returns
Source term corresponding to the user-defined exact solution, diffusivity, advection and reaction coefficients, evaluated at (x,y,z)