41 #ifndef FSIEXACTJACOBIAN_HPP 42 #define FSIEXACTJACOBIAN_HPP 44 #include <lifev/fsi/solver/FSIOperator.hpp> 70 typedef super::vector_Type vector_Type;
77 typedef super::solid_Type solid_Type;
105 const vector_Type& _res,
106 const Real _linearRelTol);
116 const vector_Type& _disp,
160 void eval (
const vector_Type& _res,
UInt status);
205 int Apply (
const Epetra_MultiVector& X, Epetra_MultiVector& Y)
const;
210 std::cout <<
"********* EJ : transpose not available\n";
216 std::cout <<
"********* EJ : inverse not available\n";
221 std::cout <<
"********* EJ : NormInf not available\n";
226 return "exactJacobian";
243 return *M_operatorDomainMap;
247 return *M_operatorRangeMap;
281 return new FSIExactJacobian();
void setOperator(FSIExactJacobian *ej)
sets the exactJacobian pointer and some contents thereof
void solveLinearSolid()
Solves the linear structure problem.
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
const Epetra_Comm & Comm() const
super::vectorPtr_Type vectorPtr_type
void registerMyProducts()
register the product for the factory
void eval(const vector_Type &_res, UInt status)
void solveLinearFluid()
Solves the linear fluid problem.
std::shared_ptr< Epetra_Comm > M_comm
FSIExactJacobian()
Empty Constructor.
void setDataFile(GetPot const &data)
initializes the GetPot data file
bool UseTranspose() const
virtual ~Epetra_ExactJacobian()
Destructor.
const Epetra_Map & OperatorDomainMap() const
std::shared_ptr< map_Type > mapPtr_Type
void bcManageVec(super::fluidBchandler_Type &, vector_Type &)
should call bcManage for a vector, but the implementation is empty
const Epetra_Map & OperatorRangeMap() const
SolverAztecOO M_linearSolver
NonLinearAitken - LifeV class for the non-linear generalized Aitken algorithm.
Epetra_ExactJacobian()
Empty Constructor.
matrixPtr_Type M_matrShapeDer
Epetra_ExactJacobian This class implements an Epetra_Operator to be passed to AztecOO.
double Real
Generic real data.
FSIModelExactJacobian - Implementation of an FSI (Operator) with Newton algorithm.
int ApplyInverse(const Epetra_MultiVector &, Epetra_MultiVector &) const
void setupFluidSolid()
setup of the fluid and solid solver classes
void setupFEspace()
sets the space discretization parameters
void evalResidual(vector_Type &_res, const vector_Type &_disp, const UInt _iter)
Evaluates the nonlinear residual of the FSI system.
int SetUseTranspose(bool)
These are the methods necessary to implement Epetra_Operator but that are not used.
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
apply the jacobian to X and returns the result in Y
NonLinearAitken< vector_Type > M_aitkFS
mapPtr_Type M_operatorRangeMap
Epetra_ExactJacobian M_epetraOper
~FSIExactJacobian()
Destructor.
fluid_Type::matrixPtr_Type matrixPtr_Type
uint32_type UInt
generic unsigned integer (used mainly for addressing)
const char * Label() const
fluid_Type::matrix_Type matrix_Type
mapPtr_Type M_operatorDomainMap
void solveJac(vector_Type &_muk, const vector_Type &_res, const Real _linearRelTol)
solves the Jacobian system