6 #ifndef PZEULERANALYSIS_H 7 #define PZEULERANALYSIS_H 49 virtual void Run(std::ostream &out,
const std::string & dxout,
int dxRes);
52 virtual void Run(std::ostream &out)
54 std::cout <<__PRETTY_FUNCTION__ <<
" should never be called!!!\n";
132 std::cout << __PRETTY_FUNCTION__ <<
" should never be called\n";
142 int RunNewton(REAL & epsilon,
int & numIter);
180 void CFLControl(REAL & lastEpsilon, REAL & epsilon, REAL & epsilon_Newton, REAL & timeStep);
void CompareRhs()
Trying difference between two times.
REAL fTimeIntEps
Stop criteria for time integration loops.
Contains TPZAnalysis class which implements the sequence of actions to perform a finite element analy...
virtual void AssembleRhs()
Assembles the right hand side vector.
virtual void Solve()
Invert the stiffness matrix.
virtual void Run(std::ostream &out, const std::string &dxout, int dxRes)
See declaration in the base class.
TPZFMatrix< STATE > fRhsLast
Vector to hold the contribution of last state to the rhs.
virtual void Assemble()
Assembles the stiffness matrix.
Implements the interface of the graphmesh to the OpenDX graphics package. Post processing.
void UpdateSolAndRhs(TPZFMatrix< STATE > &deltaSol, REAL &epsilon)
Adds deltaSol to the last solution and stores it as the current solution, preparing it to contribute ...
void WriteCMesh(const char *str)
Writes the computational mesh onto disk.
Implements an analysis procedure for computing the steady state solution of a compressible Euler flow...
int LineSearch(REAL &residual, TPZFMatrix< STATE > &sol0, TPZFMatrix< STATE > &dir)
This method will search for the solution which minimizes the residual.
int RunNewton(REAL &epsilon, int &numIter)
Implements the Newton's method.
void SetContributionTime(TPZContributeTime time)
Informs the Analysis class the time at which the current solution in the computational mesh belongs ...
int fTotalNewton
Total number of newton iterations during this run.
Declarates the TPZBlock<REAL>class which implements block matrices.
void BufferLastStateAssemble()
Buffers the assemblage of the rhs with respect to the last state.
REAL EvaluateFluxEpsilon()
Evaluates the flux part of the residual for convergence check.
int fEvolCFL
Indicates whether the CFL is to evolute or not.
REAL fNewtonEps
Stop criteria for the Newton loops.
Contains the TPZStructMatrixOR class which responsible for a interface among Matrix and Finite Elemen...
Contains TPZBlockDiagonal class which defines block diagonal matrices.
void SetEvolCFL(int EvolCFL)
Indicates whether the CFL is to evolute or not.
void SetGMResBlock(REAL tol, int numiter, int numvec)
Set to solve with GMRes algorithm for block diagonal matrix.
Computational mesh with additional data for CFD problems. Computational Mesh.
void SetFrontalSolver()
Set to solve with Frontal method.
virtual void Run(std::ostream &out)
See declaration in the base class.
void SetTimeIntCriteria(REAL epsilon, int maxIter)
Settings for the time integration method.
Implements the sequence of actions to perform a finite element analysis. Analysis.
int fNewtonMaxIter
Maxime iterations for iterative Newton loops.
void SetGMResFront(REAL tol, int numiter, int numvectors)
Set to solve with GMRes algorithm.
void SetBlockDiagonalPrecond(TPZBlockDiagonal< STATE > *blockDiag)
Informs a block diagonal to be used as preconditioning.
~TPZEulerAnalysis()
Simple destructor.
TPZFlowCompMesh * fFlowCompMesh
Stores a pointer to the computational flow mesh.
REAL ComputeTimeStep()
Defines the time step for each material in the mesh.
void UpdateHistory()
After a call to UpdateSolution, this method shifts the history in one solution, so that the current s...
Contains declaration of TPZCompMesh class which is a repository for computational elements...
virtual void Run(std::ostream &out=std::cout)
Calls the appropriate sequence of methods to build a solution or a time stepping sequence.
TPZDXGraphMesh * PrepareDXMesh(const std::string &dxout, int dxRes)
Prepares the DX graph mesh.
void CFLControl(REAL &lastEpsilon, REAL &epsilon, REAL &epsilon_Newton, REAL &timeStep)
Evaluates the CFL control based on the newest residual norm (epsilon) and last residual norm (lastEps...
Contains declaration of TPZFlowCompMesh class which is a computational mesh with additional data for ...
Contains TPZSolver class which defines a abstract class of solvers which will be used by matrix class...
int fHasFrontalPreconditioner
Indication if a frontal matrix is being used as a preconditioner.
Contains TPZStepSolver class which defines step solvers class.
virtual void Solve()
Invert the stiffness matrix.
int fTimeIntMaxIter
Maxime iterations for iterative time integration loops.
Contains the TPZConservationLaw class which implements the interface for conservation laws...
Contains the TPZDXGraphMesh class which implements the interface of the graphmesh to the OpenDX graph...
TPZBlockDiagonal< STATE > * fpBlockDiag
Preconditioner.
void SetLastState()
Sets the solution vector to be the one representing the Current State (Wn) or the last explicit solut...
TPZContributeTime
Indicates which term is put in the right hand side and tangent matrix.
void SetNewtonCriteria(REAL epsilon, int maxIter)
Settings for the Newton's method.
TPZEulerAnalysis()
Default constructor.
void SetAdvancedState()
Sets the solution vector to be the one representing the Advanced State (Wn+1) or the advanced implici...