NeoPZ
Public Member Functions | Static Public Member Functions | Static Private Member Functions | Private Attributes | List of all members
TPZMGAnalysis Class Reference

Implements multigrid analysis. TPZMGAnalysis is derived from TPZAnalysis. Analysis. More...

#include <pzmganalysis.h>

Inheritance diagram for TPZMGAnalysis:
[legend]
Collaboration diagram for TPZMGAnalysis:
[legend]

Public Member Functions

virtual ~TPZMGAnalysis ()
 Destructor. More...
 
 TPZMGAnalysis (TPZCompMesh *)
 Creates an object multigrid analysis giving a computational mesh. More...
 
void AppendMesh (TPZCompMesh *mesh)
 Append a mesh to the meshes vector. More...
 
TPZCompMeshPopMesh ()
 Pop the last mesh of the meshes vector. More...
 
virtual void Solve ()
 Uses fSolver object to apply a solution algorithm. More...
 
void ComputeError (TPZVec< REAL > &error)
 Loads the last two solutions and call the error between these two aproximations. More...
 
- Public Member Functions inherited from TPZAnalysis
 TPZAnalysis (TPZCompMesh *mesh, bool mustOptimizeBandwidth=true, std::ostream &out=std::cout)
 Create an TPZAnalysis object from one mesh pointer. More...
 
 TPZAnalysis (TPZAutoPointer< TPZCompMesh > mesh, bool mustOptimizeBandwidth=true, std::ostream &out=std::cout)
 Create an TPZAnalysis object from one mesh auto pointer object. More...
 
void SetGuiInterface (TPZAutoPointer< TPZGuiInterface > gui)
 Defines gui interface object. More...
 
TPZAutoPointer< TPZGuiInterfaceGetGuiInterface () const
 Gets gui interface object. More...
 
bool AmIKilled ()
 Returns if the process was canceled through gui interface. More...
 
virtual void SetCompMesh (TPZCompMesh *mesh, bool mustOptimizeBandwidth)
 Set the computational mesh of the analysis. More...
 
 TPZAnalysis ()
 Create an empty TPZAnalysis object. More...
 
void Write (TPZStream &buf, int withclassid) const override
 Writes this object to the TPZStream buffer. Include the classid if withclassid = true. More...
 
void Read (TPZStream &buf, void *context) override
 read objects from the stream More...
 
virtual ~TPZAnalysis (void)
 Destructor: deletes all protected dynamic allocated objects. More...
 
void CleanUp ()
 deletes all data structures More...
 
void SetRenumber (TPZAutoPointer< TPZRenumbering > renumber)
 Change the renumbering scheme. More...
 
void OptimizeBandwidth ()
 Sets the computer connection block number from the graphical connections block number otimization. More...
 
int HighestDimension ()
 Returns the dimension of the material which has the highest dimension. More...
 
void Resequence (int firstel=-1)
 Recompute the node sequence. More...
 
int ComputeNumberofLoadCases ()
 Determine the number of load cases from the material objects and return its value. More...
 
virtual void Assemble ()
 Assemble the stiffness matrix and load vector. More...
 
virtual void AssembleResidual ()
 Assemble the load vector. More...
 
TPZFMatrix< STATE > & Rhs ()
 Returns the load vector. More...
 
TPZFMatrix< STATE > & Solution ()
 Returns the solution matrix. More...
 
TPZCompMeshMesh () const
 Returns the pointer to the computational mesh. More...
 
TPZAutoPointer< TPZStructMatrixStructMatrix ()
 Returns a reference to the structural matrix. More...
 
TPZMatrixSolver< STATE > * BuildPreconditioner (EPrecond preconditioner, bool overlap)
 Define the type of preconditioner used. More...
 
void SetStep (int step)
 ste the step for post processing More...
 
void SetThreadsForError (int nthreads)
 
int GetStep ()
 
void SetTime (REAL time)
 Sets time will be used in dx files. More...
 
REAL GetTime ()
 Gets time used in dx files. More...
 
void ShowShape (const std::string &plotfile, TPZVec< int64_t > &equationindices)
 Graphic of the solution as V3DGrap visualization. More...
 
void ShowShape (const std::string &plotfile, TPZVec< int64_t > &equationindices, int matid, const TPZVec< std::string > &varname)
 Graphic of the solution as V3DGrap visualization. More...
 
void LoadShape (double dx, double dy, int64_t numelem, TPZConnect *nod)
 Make assembling and clean the load and solution vectors. More...
 
virtual void Run (std::ostream &out=std::cout)
 Calls the appropriate sequence of methods to build a solution or a time stepping sequence. More...
 
virtual void DefineGraphMesh (int dimension, const TPZVec< std::string > &scalnames, const TPZVec< std::string > &vecnames, const std::string &plotfile)
 Define GrapMesh as V3D, DX, MV or VTK depending on extension of the file. More...
 
virtual void DefineGraphMesh (int dimension, const TPZVec< std::string > &scalnames, const TPZVec< std::string > &vecnames, const TPZVec< std::string > &tensnames, const std::string &plotfile)
 Define GrapMesh as VTK with tensorial names depending on extension of the file. More...
 
virtual void DefineGraphMesh (int dimension, const std::set< int > &matids, const TPZVec< std::string > &scalnames, const TPZVec< std::string > &vecnames, const std::string &plotfile)
 Define GrapMesh as V3D, DX, MV or VTK depending on extension of the file. More...
 
virtual void DefineGraphMesh (int dimension, const std::set< int > &matids, const TPZVec< std::string > &scalnames, const TPZVec< std::string > &vecnames, const TPZVec< std::string > &tensnames, const std::string &plotfile)
 Define GrapMesh as VTK with tensorial names depending on extension of the file. More...
 
virtual void CloseGraphMesh ()
 Clean the GrapMesh vector. More...
 
TPZGraphMeshGraphMesh (int dimension)
 Defines the postprocessing parameters for the graphical grid. More...
 
virtual void PostProcess (int resolution)
 Draw solution over mesh for all dimensions. More...
 
virtual void PostProcess (int resolution, int dimension)
 Draw solution over mesh by dimension. More...
 
void IdentifyPostProcessingMatIds (int dimension, std::set< int > &matids)
 Fill mat ids with materials with provided dimension wich are not boundary conditinos or interface. More...
 
virtual void LoadSolution ()
 Load the solution into the computable grid. More...
 
virtual void LoadSolution (const TPZFMatrix< STATE > &sol)
 Load the solution into the computable mesh considering sol as Solution vector of the analysis. More...
 
TPZVec< STATE > Integrate (const std::string &varname, const std::set< int > &matids)
 Integrate the postprocessed variable name over the elements included in the set matids. More...
 
void SetExact (std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &deriv)> f)
 Sets the pointer of the exact solution function. More...
 
virtual void PostProcess (TPZVec< REAL > &loc, std::ostream &out=std::cout)
 Compute the local error over all elements and global errors in several norms and print out. More...
 
virtual void PostProcessError (TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
 Compute the local error over all elements and global errors in several norms and print out without calculating the errors of the variables for hdiv spaces. More...
 
virtual void PostProcessErrorSerial (TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
 
virtual void PostProcessErrorParallel (TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
 
void CreateListOfCompElsToComputeError (TPZAdmChunkVector< TPZCompEl *> &elvec)
 
void Print (const std::string &name, std::ostream &out)
 Print connect and solution information. More...
 
void PrintVectorByElement (std::ostream &out, TPZFMatrix< STATE > &vec, REAL tol=1.e-10)
 Print the residual vector for those elements with entry above a given tolerance. More...
 
TPZMatrixSolver< STATE > & Solver ()
 Get the solver matrix. More...
 
void AnimateRun (int64_t num_iter, int steps, TPZVec< std::string > &scalnames, TPZVec< std::string > &vecnames, const std::string &plotfile)
 Run and print the solution step by step. More...
 
void SetSolver (TPZMatrixSolver< STATE > &solver)
 Set solver matrix. More...
 
void SetStructuralMatrix (TPZAutoPointer< TPZStructMatrix > strmatrix)
 Set structural matrix as auto pointer for analysis. More...
 
void SetStructuralMatrix (TPZStructMatrix &strmatrix)
 Set structural matrix for analysis. More...
 
int ClassId () const override
 Define the class id associated with the class. More...
 
virtual void DefineElementTable (int dimension, TPZVec< int64_t > &GeoElIds, TPZVec< REAL > &points)
 Fill the computational element vector to post processing depending over geometric mesh defined. More...
 
virtual void SetTableVariableNames (TPZVec< std::string > varnames)
 Sets the names of the variables into the data structure for post processing. More...
 
virtual void PrePostProcessTable (std::ostream &out_file)
 Prepare data to print post processing and print coordinates. More...
 
virtual void PostProcessTable (std::ostream &out_file)
 Print the solution related with the computational element vector in post process. More...
 
void PostProcessTable (TPZFMatrix< REAL > &pos, std::ostream &out=std::cout)
 Compute and print the local error over all elements in data structure of post process, also compute global errors in several norms. More...
 
- Public Member Functions inherited from TPZSavable
 TPZSavable ()
 
virtual ~TPZSavable ()
 
virtual std::list< std::map< std::string, uint64_t > > VersionHistory () const
 
virtual std::pair< std::string, uint64_t > Version () const
 
virtual bool Compare (TPZSavable *copy, bool override=false)
 Compares the object for identity with the object pointed to, eventually copy the object. More...
 
virtual bool Compare (TPZSavable *copy, bool override=false) const
 Compares the object for identity with the object pointed to, eventually copy the object. More...
 
- Public Member Functions inherited from TPZRegisterClassId
template<typename T >
 TPZRegisterClassId (int(T::*)() const)
 
 TPZRegisterClassId ()=default
 

Static Public Member Functions

static TPZCompMeshUniformlyRefineMesh (TPZCompMesh *mesh, bool withP=false)
 Proceeds the uniformly h-p refinement of mesh. More...
 
static void MeshError (TPZCompMesh *fine, TPZCompMesh *coarse, TPZVec< REAL > &ervec, void(*f)(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv), TPZVec< REAL > &truervec)
 Evaluates the error between aproximation of coarse and fine meshes. More...
 
- Static Public Member Functions inherited from TPZSavable
static std::set< TPZRestoreClassBase * > & RestoreClassSet ()
 This static function guarantees that the gMap object is available when needed. More...
 
static std::map< int, TPZRestore_t > & ClassIdMap ()
 This static function guarantees that the gMap object is available when needed. More...
 
static std::pair< std::string, uint64_t > NeoPZVersion ()
 
static void Register (TPZRestoreClassBase *restore)
 
static void RegisterClassId (int classid, TPZRestore_t fun)
 
static TPZSavableCreateInstance (const int &classId)
 

Static Private Member Functions

static REAL ElementError (TPZInterpolatedElement *fine, TPZInterpolatedElement *coarse, TPZTransform<> &tr, void(*f)(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv), REAL &truerror)
 Calculates an element error based on two aproximations. More...
 

Private Attributes

TPZStack< TPZCompMesh *> fMeshes
 Contains the computational meshes of one cycle. More...
 
TPZStack< TPZFMatrix< STATE > * > fSolutions
 Contains the meshes solutions. More...
 
TPZStack< TPZMatrixSolver< STATE > * > fSolvers
 Contains the solution method applied to the mesh. More...
 
TPZStack< TPZMatrixSolver< STATE > * > fPrecondition
 Contains the preconditioner of the solution method if the solution method is a krylov method. More...
 

Additional Inherited Members

- Public Types inherited from TPZAnalysis
enum  EPrecond { EJacobi, EBlockJacobi, EElement, ENodeCentered }
 Preconditioners which can be created by objects of this class. More...
 
- Public Attributes inherited from TPZAnalysis
std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &deriv)> fExact
 Pointer to Exact solution function, it is necessary to calculating errors. More...
 
- Protected Attributes inherited from TPZAnalysis
TPZGeoMeshfGeoMesh
 Geometric Mesh. More...
 
TPZCompMeshfCompMesh
 Computational mesh. More...
 
TPZGraphMeshfGraphMesh [3]
 Graphical mesh. More...
 
TPZFMatrix< STATE > fRhs
 Load vector. More...
 
TPZFMatrix< STATE > fSolution
 Solution vector. More...
 
TPZMatrixSolver< STATE > * fSolver
 Type of solver to be applied. More...
 
TPZVec< std::string > fScalarNames [3]
 Scalar variables names - to post process. More...
 
TPZVec< std::string > fVectorNames [3]
 Vectorial variables names - to post process. More...
 
TPZVec< std::string > fTensorNames [3]
 Tensorial variables names - to post process. More...
 
int fStep
 Time step. More...
 
REAL fTime
 Time variable which is used in dx output. More...
 
int fNthreadsError
 Number of threads to be used for post-processing error. More...
 
TPZAutoPointer< TPZStructMatrixfStructMatrix
 Structural matrix. More...
 
TPZAutoPointer< TPZRenumberingfRenumber
 Renumbering scheme. More...
 
TPZAutoPointer< TPZGuiInterfacefGuiInterface
 Pointer for gui interface object. More...
 
TTablePostProcess fTable
 

Detailed Description

Implements multigrid analysis. TPZMGAnalysis is derived from TPZAnalysis. Analysis.

Definition at line 29 of file pzmganalysis.h.

Constructor & Destructor Documentation

◆ ~TPZMGAnalysis()

TPZMGAnalysis::~TPZMGAnalysis ( )
virtual

◆ TPZMGAnalysis()

TPZMGAnalysis::TPZMGAnalysis ( TPZCompMesh cmesh)

Creates an object multigrid analysis giving a computational mesh.

Definition at line 31 of file pzmganalysis.cpp.

References TPZAnalysis::fExact, fMeshes, and TPZStack< T, NumExtAlloc >::Push().

Member Function Documentation

◆ AppendMesh()

void TPZMGAnalysis::AppendMesh ( TPZCompMesh mesh)

◆ ComputeError()

void TPZMGAnalysis::ComputeError ( TPZVec< REAL > &  error)

Loads the last two solutions and call the error between these two aproximations.

Definition at line 457 of file pzmganalysis.cpp.

References fMeshes, fSolutions, MeshError(), and TPZVec< T >::NElements().

◆ ElementError()

REAL TPZMGAnalysis::ElementError ( TPZInterpolatedElement fine,
TPZInterpolatedElement coarse,
TPZTransform<> &  tr,
void(*)(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)  f,
REAL &  truerror 
)
staticprivate

◆ MeshError()

void TPZMGAnalysis::MeshError ( TPZCompMesh fine,
TPZCompMesh coarse,
TPZVec< REAL > &  ervec,
void(*)(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)  f,
TPZVec< REAL > &  truervec 
)
static

Evaluates the error between aproximation of coarse and fine meshes.

Parameters
finerefined mesh
coarsesome father mesh of fine
ervecwill return the calculated element error
fpointer function
truerveccalculates the true error between a giving a solution

Definition at line 108 of file pzmganalysis.cpp.

References TPZGeoElSide::Dimension(), TPZCompElSide::Element(), ElementError(), TPZCompMesh::ElementVec(), TPZGeoElSide::Exists(), TPZCompElSide::Exists(), test::f, TPZGeoElSide::Father2(), TPZVec< T >::Fill(), TPZCompEl::Index(), TPZCompMesh::LoadReferences(), TPZCompMesh::MaterialVec(), TPZInterpolatedElement::NConnects(), TPZCompMesh::NElements(), TPZChunkVector< T, EXP >::NElements(), numel, TPZCompEl::Reference(), TPZCompMesh::Reference(), TPZGeoElSide::Reference(), TPZGeoMesh::ResetReference(), and TPZVec< T >::Resize().

Referenced by ComputeError().

◆ PopMesh()

TPZCompMesh * TPZMGAnalysis::PopMesh ( )

◆ Solve()

void TPZMGAnalysis::Solve ( )
virtual

◆ UniformlyRefineMesh()

TPZCompMesh * TPZMGAnalysis::UniformlyRefineMesh ( TPZCompMesh mesh,
bool  withP = false 
)
static

Member Data Documentation

◆ fMeshes

TPZStack< TPZCompMesh * > TPZMGAnalysis::fMeshes
private

Contains the computational meshes of one cycle.

Definition at line 72 of file pzmganalysis.h.

Referenced by AppendMesh(), ComputeError(), PopMesh(), Solve(), TPZMGAnalysis(), and ~TPZMGAnalysis().

◆ fPrecondition

TPZStack<TPZMatrixSolver<STATE> *> TPZMGAnalysis::fPrecondition
private

Contains the preconditioner of the solution method if the solution method is a krylov method.

The preconditioner can be used as a coarse grid iteration

Definition at line 82 of file pzmganalysis.h.

Referenced by AppendMesh(), PopMesh(), Solve(), and ~TPZMGAnalysis().

◆ fSolutions

TPZStack<TPZFMatrix<STATE> *> TPZMGAnalysis::fSolutions
private

Contains the meshes solutions.

Definition at line 75 of file pzmganalysis.h.

Referenced by AppendMesh(), ComputeError(), PopMesh(), Solve(), and ~TPZMGAnalysis().

◆ fSolvers

TPZStack<TPZMatrixSolver<STATE> *> TPZMGAnalysis::fSolvers
private

Contains the solution method applied to the mesh.

Definition at line 78 of file pzmganalysis.h.

Referenced by AppendMesh(), PopMesh(), Solve(), and ~TPZMGAnalysis().


The documentation for this class was generated from the following files: