43 #include "Galeri_ConfigDefs.h"
46 #include "Epetra_MpiComm.h"
48 #include "Epetra_SerialComm.h"
50 #include "Epetra_FECrsMatrix.h"
51 #include "Epetra_FEVector.h"
53 #include "Galeri_core_Object.h"
54 #include "Galeri_core_Workspace.h"
55 #include "Galeri_grid_Loadable.h"
56 #include "Galeri_grid_Generator.h"
57 #include "Galeri_quadrature_Hex.h"
58 #include "Galeri_problem_ScalarLaplacian.h"
59 #include "Galeri_viz_MEDIT.h"
65 using namespace Galeri;
75 const double& phi_i_derx,
76 const double& phi_i_dery,
77 const double& phi_i_derz,
79 const double& phi_j_derx,
80 const double& phi_j_dery,
81 const double& phi_j_derz)
83 return(phi_i_derx * phi_j_derx +
84 phi_i_dery * phi_j_dery +
85 phi_i_derz * phi_j_derz);
95 return(-getExactSolution(
'f', x, y, z) * phi_i);
101 return(getExactSolution(
'f', x, y, z));
112 const double& y,
const double& z)
115 return(exp(x) + exp(y) + exp(z));
116 else if (what ==
'x')
118 else if (what ==
'y')
120 else if (what ==
'z')
131 int main(
int argc,
char *argv[])
134 MPI_Init(&argc,&argv);
140 Galeri::core::Workspace::setNumDimensions(3);
142 Galeri::grid::Loadable domain, boundary;
144 int numGlobalElementsX = 2 * comm.
NumProc();
145 int numGlobalElementsY = 2;
146 int numGlobalElementsZ = 2;
152 Galeri::grid::Generator::
153 getCubeWithHexs(comm, numGlobalElementsX, numGlobalElementsY, numGlobalElementsZ,
154 mx, my, mz, domain, boundary);
156 Epetra_Map matrixMap(domain.getNumGlobalVertices(), 0, comm);
162 Galeri::problem::ScalarLaplacian<Laplacian> problem(
"Hex", 1, 8);
164 problem.integrate(domain, A, RHS);
168 problem.imposeDirichletBoundaryConditions(boundary, A, RHS, LHS);
182 list.
set(
"fact: level-of-fill", 1);
189 AztecOO solver(linearProblem);
190 solver.SetAztecOption(AZ_solver, AZ_cg);
191 solver.SetPrecOperator(Prec);
192 solver.Iterate(1550, 1e-9);
195 Galeri::viz::MEDIT::write(domain,
"sol", LHS);
198 problem.computeNorms(domain, LHS);
static double getElementLHS(const double &x, const double &y, const double &z, const double &phi_i, const double &phi_i_derx, const double &phi_i_dery, const double &phi_i_derz, const double &phi_j, const double &phi_j_derx, const double &phi_j_dery, const double &phi_j_derz)
static double getExactSolution(const char &what, const double &x, const double &y, const double &z)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
int PutScalar(double ScalarConstant)
static double getBoundaryValue(const double &x, const double &y, const double &z)
virtual int Initialize()=0
Computes all it is necessary to initialize the preconditioner.
static double getElementRHS(const double &x, const double &y, const double &z, const double &phi_i)
virtual int SetParameters(Teuchos::ParameterList &List)=0
Sets all parameters for the preconditioner.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
int main(int argc, char *argv[])
Ifpack: a function class to define Ifpack preconditioners.
static Ifpack_Preconditioner * Create(EPrecType PrecType, Epetra_RowMatrix *Matrix, const int overlap=0, bool overrideSerialDefault=false)
Creates an instance of Ifpack_Preconditioner given the enum value of the preconditioner type (can not...
virtual int Compute()=0
Computes all it is necessary to apply the preconditioner.
static char getBoundaryType(const int ID, const double &x, const double &y, const double &z)
#define IFPACK_CHK_ERR(ifpack_err)