144 fRenumber = renumber;
171 virtual void Solve();
229 virtual void Run(std::ostream &out = std::cout);
244 return fGraphMesh[dimension-1];
249 virtual void PostProcess(
int resolution,
int dimension);
304 void Print(
const std::string &name , std::ostream &out );
332 static void *ThreadWork(
void *threaddata);
340 bool fStoreError =
false;
pthread_mutex_t fGetUniqueId
Mutexes (to sum error)
TPZVec< TPZCompEl * > fCompElPtr
TPZMatrixSolver< STATE > * fSolver
Type of solver to be applied.
void LoadShape(double dx, double dy, int64_t numelem, TPZConnect *nod)
Make assembling and clean the load and solution vectors.
Represents a graphical mesh used for post processing purposes. Post processing.
REAL GetTime()
Gets time used in dx files.
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
virtual void Solve()
Invert the stiffness matrix.
void SetStructuralMatrix(TPZAutoPointer< TPZStructMatrix > strmatrix)
Set structural matrix as auto pointer for analysis.
Implements a vector class which allows to use external storage provided by the user. Utility.
void Print(const std::string &name, std::ostream &out)
Print connect and solution information.
void OptimizeBandwidth()
Sets the computer connection block number from the graphical connections block number otimization...
TPZFMatrix< STATE > fSolution
Solution vector.
virtual void PostProcessTable(std::ostream &out_file)
Print the solution related with the computational element vector in post process. ...
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
int fNthreadsError
Number of threads to be used for post-processing error.
Templated vector implementation.
Contains declaration of TPZGuiInterface class.
void SetSolver(TPZMatrixSolver< STATE > &solver)
Set solver matrix.
Contains the TPZRenumbering class which defines the behavior to implementing node sequence numbering ...
TPZCompMesh * fCompMesh
Computational mesh.
Defines a class of matrix solvers. Solver.
virtual void Assemble()
Assemble the stiffness matrix and load vector.
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.
TPZCompMesh * Mesh() const
Returns the pointer to the computational mesh.
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.
virtual void PrePostProcessTable(std::ostream &out_file)
Prepare data to print post processing and print coordinates.
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.
Implements a chunk vector with free store administration. Utility.
TPZVec< std::string > fVariableNames
Contains the TPZStructMatrixOR class which responsible for a interface among Matrix and Finite Elemen...
void SetGuiInterface(TPZAutoPointer< TPZGuiInterface > gui)
Defines gui interface object.
EPrecond
Preconditioners which can be created by objects of this class.
TPZMatrixSolver< STATE > * BuildPreconditioner(EPrecond preconditioner, bool overlap)
Define the type of preconditioner used.
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.
TPZManVector< TPZManVector< REAL, 10 >, 100 > fvalues
TPZMatrixSolver< STATE > & Solver()
Get the solver matrix.
void SetExact(std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &deriv)> f)
Sets the pointer of the exact solution function.
virtual ~TPZAnalysis(void)
Destructor: deletes all protected dynamic allocated objects.
int ClassId() const override
Define the class id associated with the class.
TPZFMatrix< STATE > & Solution()
Returns the solution matrix.
Refines geometrical mesh (all the elements) num times.
void CleanUp()
deletes all data structures
TPZAutoPointer< TPZGuiInterface > fGuiInterface
Pointer for gui interface object.
TPZAdmChunkVector< TPZCompEl * > fElvec
TPZAutoPointer< TPZStructMatrix > StructMatrix()
Returns a reference to the structural matrix.
TPZAnalysis()
Create an empty TPZAnalysis object.
void SetThreadsForError(int nthreads)
void ShowShape(const std::string &plotfile, TPZVec< int64_t > &equationindices)
Graphic of the solution as V3DGrap visualization.
void SetTime(REAL time)
Sets time will be used in dx files.
TPZGraphMesh * GraphMesh(int dimension)
Defines the postprocessing parameters for the graphical grid.
Implements the sequence of actions to perform a finite element analysis. Analysis.
expr_ dx(i) *cos(expr_.val())
virtual void LoadSolution()
Load the solution into the computable grid.
TPZAutoPointer< TPZRenumbering > fRenumber
Renumbering scheme.
virtual void PostProcessErrorParallel(TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
TPZGeoMesh * fGeoMesh
Geometric Mesh.
virtual void CloseGraphMesh()
Clean the GrapMesh vector.
#define DebugStop()
Returns a message to user put a breakpoint in.
TPZVec< std::string > fScalarNames[3]
Scalar variables names - to post process.
void Read(TPZStream &buf, void *context) override
read objects from the stream
TPZVec< int64_t > fGeoElId
TPZMatrixSolver< STATE > * BuildSequenceSolver(TPZVec< int64_t > &graph, TPZVec< int64_t > &graphindex, int64_t neq, int numcolors, TPZVec< int > &colors)
Build a sequence solver based on the block graph and its colors.
Datastructure which defines postprocessing for one dimensional meshes.
Contains declaration of the TPZAutoPointer class which has Increment and Decrement actions are mutexe...
TPZVec< std::string > fVectorNames[3]
Vectorial variables names - to post process.
void IdentifyPostProcessingMatIds(int dimension, std::set< int > &matids)
Fill mat ids with materials with provided dimension wich are not boundary conditinos or interface...
pthread_mutex_t fAccessElement
Mutexes (to choose which element is next)
void Resequence(int firstel=-1)
Recompute the node sequence.
virtual void PostProcess(int resolution)
Draw solution over mesh for all dimensions.
virtual void Run(std::ostream &out=std::cout)
Calls the appropriate sequence of methods to build a solution or a time stepping sequence.
Contains TPZMatrix<TVar>class, root matrix class.
bool AmIKilled()
Returns if the process was canceled through gui interface.
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...
void SetStep(int step)
ste the step for post processing
TPZGraphMesh * fGraphMesh[3]
Graphical mesh.
TPZAutoPointer< TPZStructMatrix > fStructMatrix
Structural matrix.
virtual void AssembleResidual()
Assemble the load vector.
This class implements a geometric mesh for the pz environment. Geometry.
Implements computational mesh. Computational Mesh.
Defines the interface for saving and reading data. Persistency.
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
virtual void LoadSolution(const TPZFMatrix< STATE > &sol)
Load the solution into the computable mesh considering sol as Solution vector of the analysis...
TPZFMatrix< STATE > & Rhs()
Returns the load vector.
virtual void PostProcessErrorSerial(TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
void CreateListOfCompElsToComputeError(TPZAdmChunkVector< TPZCompEl *> &elvec)
TPZFMatrix< STATE > fRhs
Load vector.
TPZVec< std::string > fTensorNames[3]
Tensorial variables names - to post process.
int HighestDimension()
Returns the dimension of the material which has the highest dimension.
This class defines the interface to save and restore objects from TPZStream objects. Persistency.
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.
Defines the interface of a computational element. Computational Element.
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 ca...
TPZVec< REAL > fLocations
void SetRenumber(TPZAutoPointer< TPZRenumbering > renumber)
Change the renumbering scheme.
virtual void SetCompMesh(TPZCompMesh *mesh, bool mustOptimizeBandwidth)
Set the computational mesh of the analysis.
TPZAutoPointer< TPZGuiInterface > GetGuiInterface() const
Gets gui interface object.
REAL fTime
Time variable which is used in dx output.
virtual void SetTableVariableNames(TPZVec< std::string > varnames)
Sets the names of the variables into the data structure for post processing.
int ComputeNumberofLoadCases()
Determine the number of load cases from the material objects and return its value.
This class implements a reference counter mechanism to administer a dynamically allocated object...