NeoPZ
|
Implements computational mesh. Computational Mesh. More...
#include <pzcmesh.h>
Public Member Functions | |
TPZCompMesh (TPZGeoMesh *gr=0) | |
Constructor from geometrical mesh. More... | |
TPZCompMesh (TPZAutoPointer< TPZGeoMesh > &gmesh) | |
Constructor based on an autopointer to a geometric mesh. More... | |
TPZCompMesh (const TPZCompMesh ©) | |
Copy constructor. More... | |
TPZCompMesh & | operator= (const TPZCompMesh ©) |
copy the content of the mesh More... | |
virtual | ~TPZCompMesh () |
Simple Destructor. More... | |
virtual REAL | CompareMesh (int var, char *matname) |
This method will initiate the comparison between the current computational mesh and the mesh which is referenced by the geometric mesh. More... | |
void | GetRefPatches (std::set< TPZGeoEl *> &grpatch) |
Gives all patches of the mesh. More... | |
void | GetNodeToElGraph (TPZVec< int64_t > &nodtoelgraph, TPZVec< int64_t > &nodtoelgraphinde, TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindexx) |
Gives the conects graphs. More... | |
void | GetElementPatch (TPZVec< int64_t > nodtoelgraph, TPZVec< int64_t > nodtoelgraphindex, TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindex, int64_t elind, TPZStack< int64_t > &patch) |
Gives the element patch. More... | |
void | SetName (const std::string &nm) |
Set the mesh name. More... | |
void | SetDimModel (int dim) |
Set de dimension of the domain of the problem. More... | |
void | SetReference (TPZGeoMesh *gmesh) |
Sets the geometric reference mesh. More... | |
void | SetReference (TPZAutoPointer< TPZGeoMesh > &gmesh) |
Sets the geometric reference mesh. More... | |
int | Dimension () const |
Returns the dimension of the simulation. More... | |
std::string & | Name () |
Returns the mesh name. More... | |
void | CleanUp () |
Delete all the dynamically allocated data structures. More... | |
TPZMaterial * | FindMaterial (int id) |
Find the material with identity id. More... | |
void | LoadReferences () |
Map this grid in the geometric grid. More... | |
virtual void | ComputeNodElCon () |
Compute the number of elements connected to each connect object. More... | |
virtual void | ComputeNodElCon (TPZVec< int > &nelconnected) const |
Compute the number of elements connected to each connect object. More... | |
virtual void | CleanUpUnconnectedNodes () |
Delete the nodes which have no elements connected to them. More... | |
void | ComputeElGraph (TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindex) |
Computes the connectivity graph of the elements, as appropriate for the TPZRenumbering class. More... | |
void | ComputeElGraph (TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindex, std::set< int > &mat_ids) |
Computes the connectivity graph of the elements, as appropriate for the TPZRenumbering class. More... | |
virtual void | BuildCornerConnectList (std::set< int64_t > &connectindexes) const |
adds the connect indexes associated with base shape functions to the set More... | |
int | GetDefaultOrder () |
void | SetDefaultOrder (int order) |
TPZCompMesh * | Clone () const |
Clone this mesh. More... | |
void | CopyMaterials (TPZCompMesh &mesh) const |
Copies the materials of this mesh to the given mesh. More... | |
REAL | DeltaX () |
REAL | MaximumRadiusOfMesh () |
REAL | LesserEdgeOfMesh () |
TPZCompMesh * | ComputeMesh (TPZVec< int64_t > &accumlist, int64_t numaggl) |
Creates a mesh inserting into accumlist the element list to more refined mesh. More... | |
void | ComputeFillIn (int64_t resolution, TPZFMatrix< REAL > &fillin) |
This method will fill the matrix passed as parameter with a representation of the fillin of the global stiffness matrix, based on the sequence number of the connects. More... | |
int | ClassId () const override |
Returns the unique identifier for reading/writing objects to streams. More... | |
void | Write (TPZStream &buf, int withclassid) const override |
Save the element data to a stream. More... | |
void | Read (TPZStream &buf, void *context) override |
Read the element data from a stream. More... | |
Access_Private_Data | |
Routines for access to private data | |
int64_t | NConnects () const |
Number of connects allocated including free nodes. More... | |
int64_t | NIndependentConnects () |
Number of independent connect objects. More... | |
int64_t | NElements () const |
Number of computational elements allocated. More... | |
size_t | NMaterials () const |
Number of materials. More... | |
Access_Internal_Data_Structures | |
Methods for accessing the internal data structures | |
TPZAdmChunkVector< TPZCompEl * > & | ElementVec () |
Returns a reference to the element pointers vector. More... | |
TPZCompEl * | Element (int64_t iel) |
const TPZCompEl * | Element (int64_t iel) const |
const TPZAdmChunkVector< TPZCompEl * > & | ElementVec () const |
Returns a reference to the element pointers vector. More... | |
TPZAdmChunkVector< TPZConnect > & | ConnectVec () |
Return a reference to the connect pointers vector. More... | |
const TPZAdmChunkVector< TPZConnect > & | ConnectVec () const |
std::map< int,TPZMaterial *> & | MaterialVec () |
Returns a reference to the material pointers vector. More... | |
std::map< int,TPZMaterial *> | MaterialVec () const |
Returns a reference to the material pointers vector. More... | |
TPZGeoMesh * | Reference () const |
Returns a pointer to the geometrical mesh associated. More... | |
const TPZBlock< STATE > & | Block () const |
Access the block structure of the solution vector. More... | |
TPZBlock< STATE > & | Block () |
Access the block structure of the solution vector. More... | |
TPZFMatrix< STATE > & | Solution () |
Access the solution vector. More... | |
TPZFMatrix< STATE > & | SolutionN () |
Access the previous solution vector. More... | |
TPZFMatrix< STATE > & | ElementSolution () |
Access method for the element solution vectors. More... | |
Modify_Mesh_Structure | |
Methods for modify mesh structure, materials, elements, node etc | |
virtual int64_t | AllocateNewConnect (int nshape, int nstate, int order) |
Returns an index to a new connect. More... | |
virtual int64_t | AllocateNewConnect (const TPZConnect &connect) |
Returns an index to a new connect. More... | |
int | InsertMaterialObject (TPZMaterial *mat) |
Insert a material object in the datastructure. More... | |
void | InitializeBlock () |
Resequence the block object, remove unconnected connect objects and reset the dimension of the solution vector. More... | |
virtual void | ExpandSolution () |
Adapt the solution vector to new block dimensions. More... | |
Access_Solution | |
Methods for access and manipulates solution | |
void | SetElementSolution (int64_t i, TPZVec< STATE > &sol) |
Set a ith element solution, expanding the element-solution matrix if necessary. More... | |
Print | |
Print methods | |
virtual void | Print (std::ostream &out=std::cout) const |
Prints mesh data. More... | |
void | ConnectSolution (std::ostream &out) |
Print the solution by connect index. More... | |
Submesh_methods | |
Methods required by submeshes | |
virtual TPZCompMesh * | FatherMesh () const |
Get the father meshes stack. More... | |
virtual void | MakeInternal (int64_t local) |
Makes a specified connection a internal mesh connection. More... | |
virtual void | MakeAllInternal () |
Make all mesh connections internal mesh connections. Connects to an internal connection. More... | |
virtual TPZCompMesh * | RootMesh (int64_t local) |
Returns the rootmesh who have the specified connection. More... | |
virtual int64_t | TransferElementFrom (TPZCompMesh *mesh, int64_t elindex) |
Transfer one element from a specified mesh to the current submesh. More... | |
virtual int64_t | TransferElementTo (TPZCompMesh *mesh, int64_t elindex) |
Transfer one element from a submesh to another mesh. More... | |
virtual int64_t | TransferElement (TPZCompMesh *mesh, int64_t elindex) |
Transfer one element form a submesh to another mesh. More... | |
virtual int64_t | PutinSuperMesh (int64_t local, TPZCompMesh *super) |
Put an local connection in the supermesh - Supermesh is one mesh who contains the analised submesh. More... | |
virtual int64_t | GetFromSuperMesh (int64_t superind, TPZCompMesh *super) |
Get an external connection from the supermesh - Supermesh is one mesh who contains the analised submesh. More... | |
virtual TPZCompMesh * | CommonMesh (TPZCompMesh *mesh) |
Gives the commom father mesh of the specified mesh and the current submesh. 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 | |
Protected Member Functions | |
void | ModifyPermute (TPZVec< int64_t > &permute, int64_t lagrangeq, int64_t maxeq) |
Modify the permute vector swapping the lagrangeq with maxeq and shifting the intermediate equations. More... | |
Protected Attributes | |
TPZGeoMesh * | fReference |
Geometric grid to which this grid refers. More... | |
TPZAutoPointer< TPZGeoMesh > | fGMesh |
Autopointer to the geometric mesh used in case the user has passed an autopointer. More... | |
std::string | fName |
Grid name for model identification. More... | |
TPZAdmChunkVector< TPZCompEl * > | fElementVec |
List of pointers to elements. More... | |
TPZAdmChunkVector< TPZConnect > | fConnectVec |
List of pointers to nodes. More... | |
std::map< int, TPZMaterial *> | fMaterialVec |
Map of pointers to materials. More... | |
TPZBlock< STATE > | fSolutionBlock |
Block structure of the solution vector ???? More... | |
TPZFMatrix< STATE > | fSolution |
Solution vector. More... | |
TPZFMatrix< STATE > | fSolN |
Solution at previous state. More... | |
TPZBlock< STATE > | fBlock |
Block structure to right construction of the stiffness matrix and load vector. More... | |
TPZFMatrix< STATE > | fElementSolution |
Solution vectors organized by element. More... | |
int | fDimModel |
int | fDefaultOrder |
Default order for all elements of this mesh. More... | |
TPZCreateApproximationSpace | fCreate |
The object which defines the type of space being created. More... | |
int64_t | fNmeshes |
SCIENTIFIC_ROUTINES | |
int64_t | NEquations () |
This computes the number of equations associated with non-restrained nodes. More... | |
int | BandWidth () |
This method computes the bandwidth of the system of equations. More... | |
virtual void | Skyline (TPZVec< int64_t > &skyline) |
This method computes the skyline of the system of equations. More... | |
void | AssembleError (TPZFMatrix< REAL > &estimator, int errorid) |
Assemble the vector with errors estimators. More... | |
void | BuildTransferMatrix (TPZCompMesh &coarsemesh, TPZTransfer< STATE > &transfer) |
Builds the transfer matrix from the current grid to the coarse grid. More... | |
void | BuildTransferMatrixDesc (TPZCompMesh &transfermesh, TPZTransfer< STATE > &transfer) |
To discontinuous elements. More... | |
void | ProjectSolution (TPZFMatrix< STATE > &projectsol) |
void | SetNMeshes (int64_t nmeshes) |
int64_t | GetNMeshes () |
TPZCompEl * | CreateCompEl (TPZGeoEl *gel, int64_t &index) |
Create a computational element based on the geometric element. More... | |
virtual void | AutoBuild (const std::set< int > &MaterialIDs) |
Creates the computational elements, and the degree of freedom nodes. More... | |
virtual void | AutoBuild () |
Creates the computational elements, and the degree of freedom nodes. More... | |
void | AutoBuild (const TPZVec< int64_t > &gelindexes) |
build the computational elements for the geometric element indexes More... | |
virtual void | AutoBuildContDisc (const TPZVec< TPZGeoEl *> &continuous, const TPZVec< TPZGeoEl *> &discontinuous) |
Creates the computational elements, and the degree of freedom nodes. More... | |
TPZCreateApproximationSpace & | ApproxSpace () |
void | SetAllCreateFunctionsDiscontinuous () |
void | SetAllCreateFunctionsContinuous () |
void | SetAllCreateFunctionsDiscontinuousReferred () |
void | SetAllCreateFunctionsContinuousReferred () |
void | SetAllCreateFunctionsHDiv () |
void | SetAllCreateFunctionsHDivPressure () |
void | SetAllCreateFunctions (TPZCompEl &cel) |
void | SetAllCreateFunctionsMultiphysicElem () |
void | SetAllCreateFunctionsMultiphysicElemWithMem () |
void | SetAllCreateFunctionsContinuousWithMem () |
int | Consolidate () |
Will build the list of element boundary conditions build the list of connect boundary conditions. More... | |
int | Check () |
void | LoadSolution (const TPZFMatrix< STATE > &sol) |
Given the solution of the global system of equations, computes and stores the solution for the restricted nodes. More... | |
void | UpdatePreviousState (STATE mult=1.0) |
void | TransferMultiphysicsSolution () |
Transfer multiphysics mesh solution. More... | |
void | Divide (int64_t index, TPZVec< int64_t > &subindex, int interpolate=0) |
Divide the element corresponding to index. More... | |
void | Coarsen (TPZVec< int64_t > &elements, int64_t &index, bool CreateDiscontinuous=false) |
Create a computational element father for the comp. elem. into elements. More... | |
void | RemakeAllInterfaceElements () |
Deletes all interfaces and rebuild them all. More... | |
void | AdjustBoundaryElements () |
Will refine the elements associated with a boundary condition till there are no elements constrained by boundary condition elements. More... | |
void | Permute (TPZVec< int64_t > &permute) |
Permute the sequence number of the connect objects It is a permute gather operation. More... | |
void | SaddlePermute () |
Put the sequence number of the pressure connects after the seq number of the flux connects. More... | |
void | SaddlePermute2 () |
TPZVec< STATE > | Integrate (const std::string &varname, const std::set< int > &matids) |
Integrate the variable name over the mesh. More... | |
void | EvaluateError (std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)> fp, bool store_error, TPZVec< REAL > &errorSum) |
Evaluates the error given the two vectors of the analised parameters. More... | |
void | ConvertDiscontinuous2Continuous (REAL eps, int opt, int dim, TPZVec< STATE > &celJumps) |
This method compute the jump solution of interface and convert discontinuous elements with jump less than eps in continuous elements. More... | |
void | Discontinuous2Continuous (int64_t disc_index, int64_t &new_index) |
This method convert a discontinuous element with index disc_index in continuous element. More... | |
virtual void | AutoBuild (const std::set< int > *MaterialIDs) |
Creates the computational elements, and the degree of freedom nodes. More... | |
Additional Inherited Members | |
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 TPZSavable * | CreateInstance (const int &classId) |
Implements computational mesh. Computational Mesh.
The computational mesh is a repository for computational elements, nodes and material objects
The computational mesh also contains the current solution of the mesh and an elementwise solution vector
The data structure of this object is rather simple
TPZCompMesh::TPZCompMesh | ( | TPZGeoMesh * | gr = 0 | ) |
Constructor from geometrical mesh.
gr | pointer to geometrical reference mesh |
Definition at line 63 of file pzcmesh.cpp.
References TPZGeoMesh::Dimension(), fBlock, fDefaultOrder, fDimModel, fNmeshes, fReference, fSolution, fSolutionBlock, TPZCompEl::GetgOrder(), LOGPZ_DEBUG, TPZGeoMesh::Name(), TPZGeoMesh::ResetReference(), SetDimModel(), TPZBlock< TVar >::SetMatrix(), SetName(), and TPZGeoMesh::SetReference().
Referenced by Clone().
TPZCompMesh::TPZCompMesh | ( | TPZAutoPointer< TPZGeoMesh > & | gmesh | ) |
Constructor based on an autopointer to a geometric mesh.
Definition at line 99 of file pzcmesh.cpp.
References TPZGeoMesh::Dimension(), fBlock, fDefaultOrder, fDimModel, fNmeshes, fReference, fSolution, fSolutionBlock, TPZCompEl::GetgOrder(), LOGPZ_DEBUG, TPZGeoMesh::Name(), TPZGeoMesh::ResetReference(), TPZBlock< TVar >::SetMatrix(), SetName(), and TPZGeoMesh::SetReference().
TPZCompMesh::TPZCompMesh | ( | const TPZCompMesh & | copy | ) |
Copy constructor.
Update data into all the connects
Definition at line 1623 of file pzcmesh.cpp.
References TPZCompEl::Clone(), ComputeNodElCon(), CopyMaterials(), fBlock, fDefaultOrder, fDimModel, fElementVec, fName, fReference, fSolution, fSolutionBlock, LOGPZ_DEBUG, TPZChunkVector< T, EXP >::NElements(), TPZGeoMesh::ResetReference(), TPZAdmChunkVector< T, EXP >::Resize(), and TPZBlock< TVar >::SetMatrix().
|
virtual |
Simple Destructor.
Definition at line 134 of file pzcmesh.cpp.
References CleanUp(), LOGPZ_DEBUG, Print(), TPZGeoMesh::Reference(), Reference(), and TPZGeoMesh::ResetReference().
void TPZCompMesh::AdjustBoundaryElements | ( | ) |
Will refine the elements associated with a boundary condition till there are no elements constrained by boundary condition elements.
Definition at line 1423 of file pzcmesh.cpp.
References Divide(), TPZCompElSide::EqualLevelElementList(), TPZGeoElSide::Exists(), TPZGeoElSide::Father2(), pzrefine::fatherside, fElementVec, TPZCompElSide::HigherLevelElementList(), InitializeBlock(), LOGPZ_DEBUG, TPZCompEl::Material(), TPZGeoEl::NCornerNodes(), TPZGeoEl::NeighbourExists(), TPZVec< T >::NElements(), TPZChunkVector< T, EXP >::NElements(), TPZGeoEl::NSides(), porder, TPZInterpolatedElement::PreferredSideOrder(), TPZInterpolatedElement::PRefine(), TPZCompEl::Reference(), TPZManVector< T, NumExtAlloc >::Resize(), TPZGeoEl::SideDimension(), and TPZGeoEl::TypeName().
Referenced by TPZBuildMultiphysicsMesh::BuildHybridMesh(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), ConvertDiscontinuous2Continuous(), TPZAnalysisError::hp_Adaptive_Mesh_Design(), SetAllCreateFunctionsContinuousWithMem(), TPZNonLinMultGridAnalysis::UniformlyRefineMesh(), TPZBuildMultiphysicsMesh::UniformRefineCompEl(), and TPZBuildMultiphysicsMesh::UniformRefineCompMesh().
|
inlinevirtual |
Returns an index to a new connect.
nshape | Number of function shapes |
nstate | Number of variable state, depending of the material. |
order | Order of approximation |
Reimplemented in TPZSubCompMesh.
Definition at line 693 of file pzcmesh.h.
References TPZAdmChunkVector< T, EXP >::AllocateNewElement(), fBlock, fConnectVec, TPZBlock< TVar >::NBlocks(), TPZBlock< TVar >::Set(), TPZConnect::SetLagrangeMultiplier(), TPZBlock< TVar >::SetNBlocks(), TPZConnect::SetNShape(), TPZConnect::SetNState(), TPZConnect::SetOrder(), and TPZConnect::SetSequenceNumber().
Referenced by TPZSubCompMesh::AllocateNewConnect(), TPZCompElDisc::CreateMidSideConnect(), TPZAgglomerateElement::CreateMidSideConnect(), TPZInterpolatedElement::CreateMidSideConnect(), ElementSolution(), TPZSubCompMesh::MakeExternal(), TPZCreateApproximationSpace::MakeRaviartThomas(), TPZSubCompMesh::SetNumberRigidBodyModes(), hdivCurvedJCompAppMath::SetupDisconnectedHdivboud(), TPZHybridizeHDiv::SplitConnects(), TPZMHMixedHybridMeshControl::SplitFluxElementsAroundFractures(), TPZCompElHDivPressure< TSHAPE >::TPZCompElHDivPressure(), and TPZCompElHDivPressureBound< TSHAPE >::TPZCompElHDivPressureBound().
|
inlinevirtual |
Returns an index to a new connect.
connect | Connect from which attributes should be copied |
Reimplemented in TPZSubCompMesh.
Definition at line 711 of file pzcmesh.h.
References TPZAdmChunkVector< T, EXP >::AllocateNewElement(), fBlock, fConnectVec, TPZBlock< TVar >::NBlocks(), TPZConnect::NShape(), TPZConnect::NState(), TPZConnect::RemoveDepend(), TPZBlock< TVar >::Set(), TPZBlock< TVar >::SetNBlocks(), and TPZConnect::SetSequenceNumber().
|
inline |
Definition at line 498 of file pzcmesh.h.
References fCreate.
Referenced by TPZBuildSBFem::BuildComputationalMeshFromSkeleton(), TPZBuildSBFem::BuildComputationMesh(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZMultiphysicsElement::CreateInterface(), TPZMHMeshControl::CreateInternalElements(), TPZMHMixedHybridMeshControl::CreateInternalFluxElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZHybridizeHDiv::CreateMultiphysicsMesh(), CreatePressureMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CreateSkeleton(), TPZBuildSBFem::CreateVolumetricElements(), TPZBuildSBFem::CreateVolumetricElementsFromSkeleton(), TPZMHMeshControl::CriaMalhaTemporaria(), TPZCreateApproximationSpace::Hybridize(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZMHMeshControl::HybridizeSkeleton(), TPZMHMixedMeshControl::InsertPeriferalHdivMaterialObjects(), PYBIND11_MODULE(), TPZPostProcAnalysis::SetAllCreateFunctionsPostProc(), TPZReducedSpace::SetAllCreateFunctionsReducedSpace(), TPZElastoPlasticAnalysis::SetAllCreateFunctionsWithMem(), and TPZHybridizeHDiv::SplitConnects().
void TPZCompMesh::AssembleError | ( | TPZFMatrix< REAL > & | estimator, |
int | errorid | ||
) |
Assemble the vector with errors estimators.
estimator | vector where will be assembled the errors |
errorid | index for dual or wheeler estimator |
Definition at line 2113 of file pzcmesh.cpp.
References TPZCompEl::ComputeError(), TPZInterfaceElement::ComputeErrorFace(), fabs, fElementVec, TPZVec< T >::Fill(), TPZCompEl::Index(), TPZInterfaceElement::LeftElement(), NElements(), TPZVec< T >::NElements(), TPZFMatrix< TVar >::Resize(), TPZInterfaceElement::RightElement(), and TPZFMatrix< TVar >::Zero().
Referenced by SetDefaultOrder().
|
privatevirtual |
Creates the computational elements, and the degree of freedom nodes.
If MaterialIDs is passed, only element of material id in the set<int> will be created
Definition at line 308 of file pzcmesh.cpp.
References TPZCreateApproximationSpace::BuildMesh(), TPZCheckGeom::CheckIds(), fCreate, and Reference().
Referenced by TPZBuildSBFem::BuildComputationalMeshFromSkeleton(), TPZBuildSBFem::BuildComputationMesh(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), Hdiv2dPaper201504::CMeshFlux(), hdivCurvedJCompAppMath::CMeshFlux(), Hdiv3dPaper201504::CMeshFlux(), Hdiv2dPaper201504::CMeshH1(), hdivCurvedJCompAppMath::CMeshH1(), Hdiv3dPaper201504::CMeshH1(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), Input::CreateCuboSkyMatrix(), TPZHybridizeHDiv::CreateMultiphysicsMesh(), CreatePressureMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZBuildSBFem::CreateVolumetricElements(), TPZBuildSBFem::CreateVolumetricElementsFromSkeleton(), TCedricTest::GenerateCompMesh(), TPZGenSubStruct::GenerateMesh(), TPZMulticamadaOrthotropic::GenerateMesh(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZSpStructMatrix::main(), TPZSymetricSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), main(), TPZSubCompMesh::main(), PYBIND11_MODULE(), and TPZReadMeshHR::ReadMesh().
|
inlinevirtual |
Creates the computational elements, and the degree of freedom nodes.
Only element of material id in the set<int> will be created
Definition at line 469 of file pzcmesh.h.
References TPZCreateApproximationSpace::BuildMesh().
|
inlinevirtual |
Creates the computational elements, and the degree of freedom nodes.
Reimplemented in TPZFlowCompMesh, and TPZMultiphysicsCompMesh.
Definition at line 474 of file pzcmesh.h.
References TPZCreateApproximationSpace::BuildMesh(), TPZCheckGeom::CheckUniqueId(), and Reference().
Referenced by TPZFlowCompMesh::AutoBuild(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpace(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory(), and GetNMeshes().
|
inline |
build the computational elements for the geometric element indexes
Definition at line 486 of file pzcmesh.h.
References AutoBuildContDisc(), and TPZCreateApproximationSpace::BuildMesh().
|
virtual |
Creates the computational elements, and the degree of freedom nodes.
Elements created may be TPZInterpolatedElement or TPZCompElDisc.
indices contains the type of the element. Element type are given by the enumerate MCreationType.
Definition at line 326 of file pzcmesh.cpp.
References TPZCheckGeom::CheckIds(), CreateCompEl(), TPZGeoMesh::ElementVec(), fBlock, fCreate, TPZGeoEl::HasSubElement(), InitializeBlock(), TPZBlock< TVar >::NBlocks(), TPZVec< T >::NElements(), TPZChunkVector< T, EXP >::NElements(), TPZGeoEl::NumInterfaces(), TPZGeoEl::Print(), Reference(), TPZCreateApproximationSpace::SetAllCreateFunctionsContinuous(), TPZCreateApproximationSpace::SetAllCreateFunctionsDiscontinuous(), and TPZBlock< TVar >::SetNBlocks().
Referenced by AutoBuild().
int TPZCompMesh::BandWidth | ( | ) |
This method computes the bandwidth of the system of equations.
Este metodo nao e apropriado para aproximantes descontinuas
Definition at line 747 of file pzcmesh.cpp.
References TPZCompEl::BuildConnectList(), fBlock, fConnectVec, fElementVec, TPZConnect::HasDependency(), TPZConnect::IsCondensed(), NElements(), TPZVec< T >::NElements(), TPZBlock< TVar >::Position(), TPZManVector< T, NumExtAlloc >::Resize(), TPZConnect::SequenceNumber(), and TPZBlock< TVar >::Size().
Referenced by TPZSBandStructMatrix::Create(), TPZBandStructMatrix::Create(), and SetDefaultOrder().
|
inline |
Access the block structure of the solution vector.
Definition at line 213 of file pzcmesh.h.
References fBlock.
Referenced by TPZBuildMultiphysicsMesh::AddConnects(), TPZMultiphysicsCompMesh::AddConnects(), TPZInterpolatedElement::AdjustPreferredSideOrder(), TPZBuildMultiphysicsMesh::AppendConnects(), TPZPostProcAnalysis::AutoBuildDisc(), TPZCreateApproximationSpace::BuildMesh(), TPZAnalysis::BuildPreconditioner(), BuildTransferMatrix(), BuildTransferMatrixDesc(), TPZInterpolatedElement::CalcIntegral(), TPZAgglomerateElement::CalcStiff(), TPZCompElLagrange::CalcStiff(), TPZSubCompMesh::CalcStiff(), Hdiv2dPaper201504::ChangeExternalOrderConnects(), hdivCurvedJCompAppMath::ChangeExternalOrderConnects(), Hdiv3dPaper201504::ChangeExternalOrderConnects(), TPZElementMatrix::ComputeDestinationIndices(), ComputeFillIn(), TPZGenSubStruct::ComputeInternalEquationPermutation(), TPZDohrStructMatrix::ComputeInternalEquationPermutation(), TPZReducedSpace::ComputeSolution(), TPZInterpolatedElement::ComputeSolution(), TPZCompElDisc::ComputeSolution(), TPZCompElHDiv< TSHAPE >::ComputeSolutionHDiv(), TPZCompElHDivPressure< TSHAPE >::ComputeSolutionPressureHDiv(), ConnectSolution(), ConnectSolution(), TPZCompElDisc::CreateMidSideConnect(), TPZAgglomerateElement::CreateMidSideConnect(), TPZInterpolatedElement::CreateMidSideConnect(), TPZMGAnalysis::ElementError(), TPZGenSubStruct::IdentifyCornerNodes(), TPZDohrStructMatrix::IdentifyCornerNodes(), TPZGenSubStruct::IdentifyEqNumbers(), TPZDohrStructMatrix::IdentifyEqNumbers(), TPZElastoPlasticAnalysis::IdentifyEquationsToZero(), TPZSubCompMesh::InitializeEF(), TPZInterfaceElement::InitializeElementMatrix(), TPZInterpolationSpace::InterpolateSolution(), TCedricTest::InterpolationError(), TCedricTest::LoadInterpolation(), TPZCompMeshTools::LoadSolution(), TPZCondensedCompEl::LoadSolution(), TPZSBFemElementGroup::LoadSolution(), TPZCompEl::LoadSolution(), TPZSubCompMesh::LoadSolution(), TPZMultiphysicsCompMesh::LoadSolutionFromMeshes(), TPZMultiphysicsCompMesh::LoadSolutionFromMultiPhysics(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), TPZAnalysisError::MathematicaPlot(), TPZConnect::NDof(), TPZCompEl::NEquations(), TPZGenSubStruct::NInternalEq(), TPZSubCompMesh::NumInternalEquations(), TPZCompElHDiv< TSHAPE >::PRefine(), TPZConnect::Print(), TPZAgglomerateElement::ProjectSolution(), TPZGenSubStruct::ReorderInternalNodes(), TPZGenSubStruct::ReorderInternalNodes2(), TPZCompElDisc::SetDegree(), TPZCompElDisc::SetExternalShapeFunction(), TPZSubCompMesh::SetNumberRigidBodyModes(), TPZCompElHDivBound2< TSHAPE >::SetSideOrder(), TPZIntelGen< TSHAPE >::SetSideOrder(), TPZCompElHDivPressureBound< TSHAPE >::SetSideOrder(), TPZCompElHDiv< TSHAPE >::SetSideOrder(), TPZCompElDisc::SolutionX(), TPZCompElHDivPressure< TSHAPE >::TPZCompElHDivPressure(), TPZCompElHDivPressureBound< TSHAPE >::TPZCompElHDivPressureBound(), TPZBuildMultiphysicsMesh::TransferFromMeshes(), TPZBuildMultiphysicsMesh::TransferFromMultiPhysics(), and TPZMultiphysicsElement::TransferMultiphysicsElementSolution().
|
inline |
|
virtual |
adds the connect indexes associated with base shape functions to the set
Reimplemented in TPZSubCompMesh.
Definition at line 2782 of file pzcmesh.cpp.
References TPZCompEl::BuildCornerConnectList(), ElementVec(), and NElements().
Referenced by ElementSolution(), and TPZDohrStructMatrix::IdentifyCornerNodes().
void TPZCompMesh::BuildTransferMatrix | ( | TPZCompMesh & | coarsemesh, |
TPZTransfer< STATE > & | transfer | ||
) |
Builds the transfer matrix from the current grid to the coarse grid.
coarsemesh | grid for where the matrix will be transfered |
transfer | transfer matrix between the current mesh and the coarse mesh |
Definition at line 883 of file pzcmesh.cpp.
References Block(), TPZInterpolationSpace::BuildTransferMatrix(), TPZGeoEl::BuildTransform2(), TPZMaterial::Dimension(), TPZCompEl::Dimension(), TPZGeoEl::Father(), fElementVec, fMaterialVec, LoadReferences(), TPZCompEl::Mesh(), NElements(), NIndependentConnects(), NMaterials(), TPZGeoEl::NSides(), TPZMaterial::NStateVariables(), TPZInterpolationSpace::Print(), PZError, TPZCompEl::Reference(), Reference(), TPZGeoEl::Reference(), TPZGeoMesh::ResetReference(), and TPZTransfer< TVar >::SetBlocks().
Referenced by TPZMGAnalysis::AppendMesh(), and SetDefaultOrder().
void TPZCompMesh::BuildTransferMatrixDesc | ( | TPZCompMesh & | transfermesh, |
TPZTransfer< STATE > & | transfer | ||
) |
To discontinuous elements.
Definition at line 941 of file pzcmesh.cpp.
References Block(), TPZCompElDisc::BuildTransferMatrix(), DebugStop, TPZMaterial::Dimension(), TPZCompEl::Dimension(), EAgglomerate, EDiscontinuous, ElementVec(), fMaterialVec, LoadReferences(), TPZChunkVector< T, EXP >::NElements(), NIndependentConnects(), TPZAgglomerateElement::NIndexes(), NMaterials(), TPZMaterial::NStateVariables(), PZError, Reference(), TPZGeoMesh::ResetReference(), TPZTransfer< TVar >::SetBlocks(), TPZAgglomerateElement::SubElement(), TPZCompEl::Type(), and TPZCompElDisc::Type().
Referenced by SetDefaultOrder().
int TPZCompMesh::Check | ( | ) |
??
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
overridevirtual |
Returns the unique identifier for reading/writing objects to streams.
returns the unique identifier for reading/writing objects to streams
Implements TPZSavable.
Reimplemented in TPZSubCompMesh, TPZFlowCompMesh, TPZCompMeshReferred, and TPZAgglomerateMesh.
Definition at line 1969 of file pzcmesh.cpp.
References Hash().
Referenced by TPZCompMeshReferred::ClassId(), TPZFlowCompMesh::ClassId(), TPZSubCompMesh::ClassId(), and SetAllCreateFunctionsContinuousWithMem().
void TPZCompMesh::CleanUp | ( | ) |
Delete all the dynamically allocated data structures.
Definition at line 161 of file pzcmesh.cpp.
References TPZAdmChunkVector< T, EXP >::CompactDataStructure(), cond, fBlock, fConnectVec, fElementVec, fMaterialVec, fSolution, fSolutionBlock, LoadReferences(), NElements(), NMaterials(), TPZAdmChunkVector< T, EXP >::NOW, TPZFMatrix< TVar >::Redim(), Reference(), TPZGeoMesh::ResetReference(), TPZAdmChunkVector< T, EXP >::Resize(), TPZBlock< TVar >::SetNBlocks(), TPZCondensedCompEl::Unwrap(), and TPZElementGroup::Unwrap().
Referenced by Name(), operator=(), ~TPZCompMesh(), and TPZMHMixedMeshControl::~TPZMHMixedMeshControl().
|
virtual |
Delete the nodes which have no elements connected to them.
Definition at line 498 of file pzcmesh.cpp.
References ComputeNodElCon(), DebugStop, fBlock, fConnectVec, TPZConnect::HasDependency(), TPZConnect::IsCondensed(), LOGPZ_DEBUG, TPZBlock< TVar >::NBlocks(), NConnects(), TPZConnect::NElConnected(), TPZChunkVector< T, EXP >::NElements(), Permute(), PZError, TPZConnect::Reset(), TPZConnect::SequenceNumber(), TPZBlock< TVar >::Set(), TPZAdmChunkVector< T, EXP >::SetFree(), and TPZBlock< TVar >::SetNBlocks().
Referenced by TPZSubCompMesh::Assemble(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), BuildElementGroups(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZSubCompMesh::CalcStiff(), Hdiv2dPaper201504::ChangeExternalOrderConnects(), hdivCurvedJCompAppMath::ChangeExternalOrderConnects(), Hdiv3dPaper201504::ChangeExternalOrderConnects(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), ConvertDiscontinuous2Continuous(), TPZCompMeshTools::CreatedCondensedElements(), TPZMHMixedMeshControl::CreateHDivPressureMHMMesh(), ElementSolution(), TPZMHMixedMeshControl::GroupandCondenseElements(), TPZMHMixedHybridMeshControl::GroupandCondenseElements(), TPZMHMixedMeshChannelControl::HideTheElements(), TPZMHMixedMeshControl::HideTheElements(), TPZCreateApproximationSpace::Hybridize(), InitializeBlock(), Hdiv2dPaper201504::PrintErrors(), hdivCurvedJCompAppMath::PrintErrors(), Hdiv3dPaper201504::PrintErrors(), ResequenceByGeometry(), Hdiv2dPaper201504::Run(), hdivCurvedJCompAppMath::Run(), Hdiv3dPaper201504::Run(), TPZDohrStructMatrix::SubStructure(), SubStructure(), TPZGenSubStruct::SubStructure(), TPZMHMeshControl::SubStructure(), TPZCreateApproximationSpace::UndoMakeRaviartThomas(), TPZBuildMultiphysicsMesh::UniformRefineCompEl(), and TPZBuildMultiphysicsMesh::UniformRefineCompMesh().
TPZCompMesh * TPZCompMesh::Clone | ( | ) | const |
Clone this mesh.
Definition at line 1720 of file pzcmesh.cpp.
References TPZCompMesh().
Referenced by SetAllCreateFunctionsContinuousWithMem().
void TPZCompMesh::Coarsen | ( | TPZVec< int64_t > & | elements, |
int64_t & | index, | ||
bool | CreateDiscontinuous = false |
||
) |
Create a computational element father for the comp. elem. into elements.
elements | indices of subelements to be agrouped |
index | of new coarse element |
CreateDiscontinuous | indicates a TPZInterpolatedElement must be created. True indicates a TPZCompElDisc |
Definition at line 1164 of file pzcmesh.cpp.
References CreateCompEl(), TPZInterpolatedElement::EDelete, ElementVec(), TPZGeoEl::Father(), fCreate, fElementVec, GetDefaultOrder(), TPZVec< T >::NElements(), TPZGeoEl::NSubElements(), TPZCompEl::Reference(), TPZInterpolationSpace::RemoveInterfaces(), TPZInterpolatedElement::RemoveSideRestraintsII(), TPZGeoEl::ResetReference(), TPZCreateApproximationSpace::SetAllCreateFunctionsContinuous(), TPZCreateApproximationSpace::SetAllCreateFunctionsDiscontinuous(), and TPZCompElDisc::SetDegree().
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
virtual |
Gives the commom father mesh of the specified mesh and the current submesh.
mesh | pointer to other mesh whosw want to know the commom father mesh |
Reimplemented in TPZSubCompMesh.
Definition at line 2794 of file pzcmesh.cpp.
References FatherMesh(), and TPZStack< T, NumExtAlloc >::Push().
Referenced by TPZInterfaceElement::ConnectIndex(), TPZInterfaceElement::GetConnects(), and TransferElement().
|
virtual |
This method will initiate the comparison between the current computational mesh and the mesh which is referenced by the geometric mesh.
var | state variable number |
matname | pointer to material name |
Definition at line 1522 of file pzcmesh.cpp.
References TPZCompEl::CompareElement(), error(), fElementVec, and TPZChunkVector< T, EXP >::NElements().
Referenced by TPZSubCompMesh::CompareElement().
void TPZCompMesh::ComputeElGraph | ( | TPZStack< int64_t > & | elgraph, |
TPZVec< int64_t > & | elgraphindex | ||
) |
Computes the connectivity graph of the elements, as appropriate for the TPZRenumbering class.
elgraph | stack of elements to create the grapho???? |
elgraphindex | graphos indexes vector |
Definition at line 1091 of file pzcmesh.cpp.
References fMaterialVec.
Referenced by TPZAnalysis::BuildPreconditioner(), ComputeFillIn(), TPZDohrStructMatrix::Create(), ElementSolution(), GetNodeToElGraph(), TPZGenSubStruct::IdentifyCornerNodes(), TPZDohrStructMatrix::IdentifyCornerNodes(), TPZAnalysis::OptimizeBandwidth(), and TPZDohrStructMatrix::SubStructure().
void TPZCompMesh::ComputeElGraph | ( | TPZStack< int64_t > & | elgraph, |
TPZVec< int64_t > & | elgraphindex, | ||
std::set< int > & | mat_ids | ||
) |
Computes the connectivity graph of the elements, as appropriate for the TPZRenumbering class.
elgraph | stack of elements to create the grapho???? |
elgraphindex | graphos indexes vector |
Apply the filter by material identifiers
Definition at line 1100 of file pzcmesh.cpp.
References TPZCompEl::BuildConnectList(), DebugStop, fConnectVec, fElementVec, TPZConnect::HasDependency(), TPZCompEl::HasMaterial(), TPZConnect::IsCondensed(), NElements(), TPZVec< T >::NElements(), NIndependentConnects(), TPZStack< T, NumExtAlloc >::Push(), TPZManVector< T, NumExtAlloc >::Resize(), TPZVec< T >::Resize(), and TPZConnect::SequenceNumber().
void TPZCompMesh::ComputeFillIn | ( | int64_t | resolution, |
TPZFMatrix< REAL > & | fillin | ||
) |
This method will fill the matrix passed as parameter with a representation of the fillin of the global stiffness matrix, based on the sequence number of the connects.
resolution | Number of rows and columns of the matrix |
fillin | Matrix which is mapped onto the global system of equations and represents the fillin be assigning a value between 0. and 1. in each element |
This method will fill the matrix passed as parameter with a representation of the fillin of the global stiffness matrix, based on the sequence number of the connects
resolution | Number of rows and columns of the matrix |
fillin | Matrix which is mapped onto the global system of equations and represents the fillin be assigning a value between 0. and 1. in each element |
Definition at line 1869 of file pzcmesh.cpp.
References Block(), ComputeElGraph(), ComputeNodElCon(), TPZRenumbering::ConvertGraph(), fConnectVec, fElementVec, TPZConnect::HasDependency(), TPZConnect::IsCondensed(), TPZVec< T >::NElements(), TPZChunkVector< T, EXP >::NElements(), NEquations(), TPZBlock< TVar >::Position(), TPZFMatrix< TVar >::Redim(), TPZConnect::SequenceNumber(), and TPZBlock< TVar >::Size().
Referenced by AssembleMatrices(), TPZDohrStructMatrix::Create(), SetAllCreateFunctionsContinuousWithMem(), TPZSubCompMesh::SetAnalysisFrontal(), and TPZSubCompMesh::SetAnalysisSkyline().
TPZCompMesh* TPZCompMesh::ComputeMesh | ( | TPZVec< int64_t > & | accumlist, |
int64_t | numaggl | ||
) |
Creates a mesh inserting into accumlist the element list to more refined mesh.
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
virtual |
Compute the number of elements connected to each connect object.
Reimplemented in TPZSubCompMesh.
Definition at line 668 of file pzcmesh.cpp.
References TPZCompEl::BuildConnectList(), fConnectVec, fElementVec, TPZConnect::IncrementElConnected(), NConnects(), NElements(), TPZVec< T >::NElements(), TPZConnect::ResetElConnected(), TPZManVector< T, NumExtAlloc >::Resize(), and TPZConnect::SequenceNumber().
Referenced by TPZFrontStructMatrix< front >::AdjustSequenceNumbering(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedHybridMeshControl::BuildComputationalMesh(), TPZMHMixedMeshControl::BuildComputationalMesh(), BuildElementGroups(), CleanUpUnconnectedNodes(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), ComputeFillIn(), TPZSubCompMesh::ComputeNodElCon(), TPZDohrStructMatrix::Create(), TCedricTest::CreateCondensedElements(), TPZMHMixedMeshControl::CreateHDivMHMMesh(), TPZMHMixedMeshControl::CreateHDivPressureMHMMesh(), ElementSolution(), TPZHybridizeHDiv::GroupandCondenseElements(), TPZMHMixedHybridMeshControl::GroupandCondenseElements(), TPZCompMeshTools::GroupElements(), TPZMHMixedHybridMeshControl::GroupElements(), TPZMHMixedMeshChannelControl::HideTheElements(), TPZMHMixedMeshControl::HideTheElements(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZGenSubStruct::IdentifyCornerNodes(), TPZGenSubStruct::InitializeDohr(), TPZGenSubStruct::InitializeDohrCondense(), TPZSubCompMesh::PotentialInternal(), TPZCompMeshTools::PutinSubmeshes(), TPZDohrStructMatrix::SubStructure(), SubStructure(), TPZGenSubStruct::SubStructure(), TPZMHMeshControl::SubStructure(), and TPZCompMesh().
|
virtual |
Compute the number of elements connected to each connect object.
Reimplemented in TPZSubCompMesh.
Definition at line 697 of file pzcmesh.cpp.
References TPZCompEl::BuildConnectList(), fElementVec, TPZVec< T >::Fill(), NConnects(), NElements(), TPZVec< T >::NElements(), TPZManVector< T, NumExtAlloc >::Resize(), and TPZVec< T >::Resize().
void TPZCompMesh::ConnectSolution | ( | std::ostream & | out | ) |
Print the solution by connect index.
out | indicates the device where the data will be printed |
Definition at line 1373 of file pzcmesh.cpp.
References Block(), ConnectVec(), NConnects(), TPZConnect::NElConnected(), TPZBlock< TVar >::Position(), TPZConnect::SequenceNumber(), TPZBlock< TVar >::Size(), and Solution().
Referenced by ElementSolution(), and TPZAnalysis::Print().
|
inline |
Return a reference to the connect pointers vector.
Definition at line 198 of file pzcmesh.h.
References fConnectVec.
Referenced by TPZBuildMultiphysicsMesh::AddConnects(), TPZMultiphysicsCompMesh::AddConnects(), TPZConnect::AddToList(), TPZCompMeshTools::AdjustFluxPolynomialOrders(), TPZFrontStructMatrix< front >::AdjustSequenceNumbering(), TPZBuildMultiphysicsMesh::AppendConnects(), TPZElementMatrix::ApplyConstraints(), AssembleMatrices(), TPZPostProcAnalysis::AutoBuildDisc(), TPZConnect::BuildConnectList(), TPZCompEl::BuildConnectList(), TPZElementMatrix::BuildDependencyOrder(), TPZConnect::BuildDependencyOrder(), TPZInterpolationSpace::BuildTransferMatrix(), TPZInterpolatedElement::BuildTransferMatrix(), TPZCompElDisc::BuildTransferMatrix(), TPZSubCompMesh::CalcStiff(), TPZConnect::CheckDependency(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZElementMatrix::ComputeDestinationIndices(), TPZGenSubStruct::ComputeInternalEquationPermutation(), TPZDohrStructMatrix::ComputeInternalEquationPermutation(), TPZSubCompMesh::ComputePermutationInternalFirst(), TPZCompEl::Connect(), ConnectSolution(), ConnectSolution(), TPZMHMeshControl::CreateInternalElements(), TPZCompElDisc::CreateMidSideConnect(), TPZAgglomerateElement::CreateMidSideConnect(), TPZInterpolatedElement::CreateMidSideConnect(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZConnect::DependencyDepth(), TPZDXGraphMesh::DrawSolution(), TPZInterpolationSpace::ExpandShapeFunctions(), GatherDependency(), TPZInterfaceElement::GetConnects(), TPZSubCompMesh::GetFromSuperMesh(), TPZElementMatrix::HasDependency(), TPZGenSubStruct::IdentifyCornerNodes(), TPZDohrStructMatrix::IdentifyCornerNodes(), TPZGenSubStruct::IdentifyEqNumbers(), TPZDohrStructMatrix::IdentifyEqNumbers(), TPZMultiphysicsInterfaceElement::IncrementElConnected(), TPZInterfaceElement::InitializeElementMatrix(), TCedricTest::InterpolationError(), TPZCompEl::LoadSolution(), TPZSubCompMesh::LoadSolution(), TPZMultiphysicsCompMesh::LoadSolutionFromMeshes(), TPZMultiphysicsCompMesh::LoadSolutionFromMultiPhysics(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), TPZSubCompMesh::MakeAllInternal(), TPZSubCompMesh::MakeExternal(), TPZSubCompMesh::MakeInternalFast(), TPZAnalysisError::MathematicaPlot(), TPZGenSubStruct::NInternalEq(), TPZFrontStructMatrix< front >::OrderElement(), TPZStructMatrixGCTP::OrderElement(), TPZStructMatrixGC::OrderElement(), TPZStructMatrixOT::OrderElement(), TPZSubCompMesh::PotentialInternal(), TPZCompEl::PressureConnectIndex(), TPZElementMatrix::Print(), TPZAnalysis::Print(), TPZAnalysis::PrintVectorByElement(), TPZGenSubStruct::ReorderInternalNodes(), TPZGenSubStruct::ReorderInternalNodes2(), SaddlePermute(), TPZCompElDisc::SetDegree(), TPZConnect::SetDependenceOrder(), TPZCompElDisc::SetExternalShapeFunction(), TPZSubCompMesh::SetNumberRigidBodyModes(), TPZCompMeshTools::SetPressureOrders(), TPZInterpolatedElement::SideConnect(), Skyline(), TPZHybridizeHDiv::SplitConnects(), TPZMHMixedHybridMeshControl::SplitFluxElementsAroundFractures(), TPZCompElHDiv< TSHAPE >::TPZCompElHDiv(), TPZCompElHDivBound2< TSHAPE >::TPZCompElHDivBound2(), TPZCompElHDivPressure< TSHAPE >::TPZCompElHDivPressure(), TPZCompElHDivPressureBound< TSHAPE >::TPZCompElHDivPressureBound(), TPZCompElLagrange::TPZCompElLagrange(), TPZIntelGen< TSHAPE >::TPZIntelGen(), TPZSubCompMesh::TransferDependencies(), TPZBuildMultiphysicsMesh::TransferFromMeshes(), TPZBuildMultiphysicsMesh::TransferFromMultiPhysics(), TPZMHMeshControl::TransferToMultiphysics(), TPZBuildMultiphysicsMesh::UniformRefineCompEl(), TPZBuildMultiphysicsMesh::UniformRefineCompMesh(), TPZIntelGen< TSHAPE >::~TPZIntelGen(), and TPZMHMixedMeshControl::~TPZMHMixedMeshControl().
|
inline |
Definition at line 200 of file pzcmesh.h.
References fConnectVec.
int TPZCompMesh::Consolidate | ( | ) |
Will build the list of element boundary conditions build the list of connect boundary conditions.
Put material pointers into the elements. Check on the number of dof of the connects
Referenced by SetAllCreateFunctionsContinuousWithMem().
void TPZCompMesh::ConvertDiscontinuous2Continuous | ( | REAL | eps, |
int | opt, | ||
int | dim, | ||
TPZVec< STATE > & | celJumps | ||
) |
This method compute the jump solution of interface and convert discontinuous elements with jump less than eps in continuous elements.
It may be compared the following values to eps:
If then
If then
eps | Tolerancy of the jump to cancel and to convert discontinuous element in continuous |
opt | Option by type of norm ( norm or norm). |
dim | Dimension of the working discontinuous elements |
celJumps | Vector to store the diference between the values from right and left elements connected on the interface |
Definition at line 2052 of file pzcmesh.cpp.
References abs(), AdjustBoundaryElements(), CleanUpUnconnectedNodes(), TPZMatrix< TVar >::Cols(), DebugStop, TPZGeoEl::Dimension(), Discontinuous2Continuous(), ElementVec(), TPZInterfaceElement::EvaluateInterfaceJump(), TPZVec< T >::Fill(), TPZCompEl::Index(), TPZInterfaceElement::LeftElement(), LOGPZ_FATAL, NElements(), TPZVec< T >::NElements(), TPZCompEl::Reference(), TPZVec< T >::Resize(), TPZInterfaceElement::RightElement(), Solution(), and sqrt.
Referenced by SetAllCreateFunctionsContinuousWithMem().
void TPZCompMesh::CopyMaterials | ( | TPZCompMesh & | mesh | ) | const |
Copies the materials of this mesh to the given mesh.
Definition at line 1797 of file pzcmesh.cpp.
References fMaterialVec.
Referenced by TPZAgglomerateElement::CreateAgglomerateMesh(), TPZHybridizeHDiv::CreateMultiphysicsMesh(), operator=(), SetAllCreateFunctionsContinuousWithMem(), TPZCompMesh(), and TPZMGAnalysis::UniformlyRefineMesh().
Create a computational element based on the geometric element.
Definition at line 462 of file pzcmesh.h.
References TPZCreateApproximationSpace::CreateCompEl().
Referenced by TPZMHMeshControl::AddBoundaryElements(), TPZBuildMultiphysicsMesh::AddWrap(), AutoBuildContDisc(), TPZPostProcAnalysis::AutoBuildDisc(), Coarsen(), TPZMHMixedHybridMeshControl::CreateAxialFluxElement(), TPZGeoEl::CreateBCCompEl(), TPZMHMixedHybridMeshControl::CreateHDivWrappers(), TPZMHMeshControl::CreateInternalElements(), TPZMHMixedHybridMeshControl::CreateInternalFluxElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZMHMixedHybridMeshControl::CreateLowerDimensionPressureElements(), TPZMHMixedHybridMeshControl::CreatePressureInterfaces(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CreateSkeleton(), Discontinuous2Continuous(), TPZCompElDisc::Divide(), TPZInterpolatedElement::Divide(), TPZMHMixedMeshControl::HybridizeSkeleton(), TPZMHMeshControl::HybridizeSkeleton(), TPZMHMeshControl::TransferToMultiphysics(), TPZMGAnalysis::UniformlyRefineMesh(), and TPZNonLinMultGridAnalysis::UniformlyRefineMesh().
REAL TPZCompMesh::DeltaX | ( | ) |
Definition at line 1810 of file pzcmesh.cpp.
References TPZGeoNode::Coord(), dist(), TPZGeoEl::Distance(), EInterface, ElementVec(), TPZMaterial::Id(), TPZCompEl::Material(), TPZChunkVector< T, EXP >::NElements(), TPZGeoEl::NodePtr(), TPZCompEl::Reference(), and TPZCompEl::Type().
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
inline |
Returns the dimension of the simulation.
Definition at line 148 of file pzcmesh.h.
References fDimModel.
Referenced by TPZBuildMultiphysicsMesh::AddWrap(), TPZCompMeshTools::AdjustFluxPolynomialOrders(), TPZNonLinMultGridAnalysis::AgglomerateMesh(), TPZHybridizeHDiv::AssociateElements(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedHybridMeshControl::BuildComputationalMesh(), TPZMHMixedMeshControl::BuildComputationalMesh(), TPZBuildSBFem::BuildComputationalMeshFromSkeleton(), TPZBuildSBFem::BuildComputationMesh(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZSBFemElementGroup::CalcStiff(), Hdiv2dPaper201504::ChangeExternalOrderConnects(), hdivCurvedJCompAppMath::ChangeExternalOrderConnects(), Hdiv3dPaper201504::ChangeExternalOrderConnects(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZCompMeshTools::ComputeDifferenceNorm(), TPZAgglomerateElement::CreateAgglomerateMesh(), TPZBuildSBFem::CreateElementGroups(), TPZMHMixedMeshControl::CreateHDivMHMMesh(), TPZHybridizeHDiv::CreateInterfaceElements(), TPZMultiphysicsElement::CreateInterfaces(), TPZCompElDisc::CreateMidSideConnect(), CreateNoElement(), CreatePressureMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CreateSkeleton(), Hdiv2dPaper201504::ErrorH1(), hdivCurvedJCompAppMath::ErrorH1(), Hdiv3dPaper201504::ErrorH1(), Hdiv2dPaper201504::ErrorPrimalDual(), hdivCurvedJCompAppMath::ErrorPrimalDual(), Hdiv3dPaper201504::ErrorPrimalDual(), TPZMultiphysicsCompEl< TGeometry >::EvaluateError(), TPZInterpolationSpace::EvaluateError(), TPZElementGroup::EvaluateError(), TPZSBFemVolume::EvaluateError(), TPZCompMeshTools::GroupElements(), TPZMHMixedHybridMeshControl::GroupElements(), TPZMHMeshControl::HybridizeSkeleton(), TPZHybridizeHDiv::InsertPeriferalMaterialObjects(), TPZInterpolationSpace::IntegrateSolution(), TPZAgglomerateElement::ListOfGroupings(), main(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), TPZAnalysis::PostProcess(), TPZAnalysis::PostProcessErrorSerial(), Print(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), TCedricTest::Run(), SetAllCreateFunctionsHDiv(), SetAllCreateFunctionsHDivPressure(), TPZGradientReconstruction::TPZGradientData::SetCel(), TPZCompMeshTools::SetPressureOrders(), TPZBuildMultiphysicsMesh::ShowShape(), TPZAnalysis::ShowShape(), TPZHybridizeHDiv::SplitConnects(), TPZDohrStructMatrix::SubStructure(), TPZCompElHDivPressureBound< TSHAPE >::TPZCompElHDivPressureBound(), TPZMeshSolution::TPZMeshSolution(), TPZSubCompMesh::TPZSubCompMesh(), TPZNonLinMultGridAnalysis::TwoGridAlgorithm(), and TPZHybridizeHDiv::VerifySolutionConsistency().
void TPZCompMesh::Discontinuous2Continuous | ( | int64_t | disc_index, |
int64_t & | new_index | ||
) |
This method convert a discontinuous element with index disc_index in continuous element.
disc_index | Index of the discontinuous element to be converted |
new_index | Returns the index of the new continuous element created |
Definition at line 1234 of file pzcmesh.cpp.
References CreateCompEl(), TPZInterpolationSpace::CreateInterfaces(), ExpandSolution(), fCreate, fElementVec, TPZInterpolationSpace::InterpolateSolution(), LOGPZ_FATAL, TPZCompEl::Reference(), TPZGeoEl::Reference(), TPZInterpolationSpace::RemoveInterfaces(), TPZGeoEl::ResetReference(), and TPZCreateApproximationSpace::SetAllCreateFunctionsContinuous().
Referenced by ConvertDiscontinuous2Continuous(), and SetAllCreateFunctionsContinuousWithMem().
void TPZCompMesh::Divide | ( | int64_t | index, |
TPZVec< int64_t > & | subindex, | ||
int | interpolate = 0 |
||
) |
Divide the element corresponding to index.
index | element to divide index |
subindex | divided elements indexes vector |
interpolate | ???? |
Definition at line 1146 of file pzcmesh.cpp.
References TPZCompEl::Divide(), fElementVec, TPZCompEl::Index(), PZError, and TPZVec< T >::Resize().
Referenced by AdjustBoundaryElements(), and SetAllCreateFunctionsContinuousWithMem().
|
inline |
Definition at line 185 of file pzcmesh.h.
Referenced by TPZSBFemElementGroup::AddElement(), TPZBuildMultiphysicsMesh::AddElements(), TPZMultiphysicsCompMesh::AddElements(), TPZCompMeshTools::AddHDivPyramidRestraints(), TPZCompMeshTools::AdjustFluxPolynomialOrders(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedMeshControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), TPZSBFemVolume::BuildCornerConnectList(), TPZMHMixedMeshControl::BuildMultiPhysicsMesh(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory(), TPZMHMixedHybridMeshControl::CheckMeshConsistency(), TPZMultiphysicsCompMesh::CleanElementsConnects(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), TPZBuildMultiphysicsMesh::ComputeAtomicIndexes(), TPZCompMeshTools::ComputeDifferenceNorm(), TPZSBFemVolume::ComputeKMatrices(), TPZSBFemVolume::ComputeSolution(), TPZMHMixedHybridMeshControl::CreateAxialFluxElement(), CreateCompMesh(), TPZCompMeshTools::CreatedCondensedElements(), TPZBuildSBFem::CreateElementGroups(), TPZMHMixedHybridMeshControl::CreateHDivWrappers(), TPZHybridizeHDiv::CreateInterfaceElements(), TPZMHMixedHybridMeshControl::CreateInternalAxialFluxes(), TPZMHMeshControl::CreateInternalElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedHybridMeshControl::CreateLowerDimensionPressureElements(), TPZMHMixedHybridMeshControl::CreateMultiPhysicsInterfaceElements(), TPZMHMixedMeshControl::CreateMultiPhysicsInterfaceElements(), TPZMHMixedHybridMeshControl::CreatePressureInterfaces(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedHybridMeshControl::CreateSkeletonAxialFluxes(), TPZMHMixedMeshControl::DeletePressureElements(), TPZStructMatrixGC::ElementColoring(), TPZStructMatrixOT::ElementColoring(), TPZCompMeshTools::ExpandHDivPyramidRestraints(), TPZHybridizeHDiv::GroupandCondenseElements(), TPZMHMixedMeshControl::GroupandCondenseElements(), TPZMHMixedHybridMeshControl::GroupandCondenseElements(), TPZCompMeshTools::GroupElements(), TPZMHMixedHybridMeshControl::GroupElements(), TPZMHMixedMeshChannelControl::HideTheElements(), TPZMHMixedMeshControl::HideTheElements(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZMHMixedMeshControl::HybridizeSkeleton(), TPZMHMeshControl::HybridizeSkeleton(), TPZAnalysis::Integrate(), Integrate(), TPZCompMeshTools::LoadSolution(), TPZStructMatrixGC::OrderElement(), TPZSBFemVolume::PRefine(), Hdiv2dPaper201504::PrintErrors(), hdivCurvedJCompAppMath::PrintErrors(), Hdiv3dPaper201504::PrintErrors(), TPZMHMixedHybridMeshControl::PrintFriendly(), TPZAnalysis::PrintVectorByElement(), TPZCompMeshTools::PutinSubmeshes(), Hdiv2dPaper201504::Run(), hdivCurvedJCompAppMath::Run(), Hdiv3dPaper201504::Run(), TPZSBFemVolume::SetElementGroupIndex(), TPZSBFemVolume::SetPreferredOrder(), TPZCompMeshTools::SetPressureOrders(), TPZSBFemVolume::SetSkeleton(), hdivCurvedJCompAppMath::SetupDisconnectedHdivboud(), TPZSBFemVolume::Shape(), TPZMHMixedHybridMeshControl::SplitFluxElementsAroundFractures(), TPZMHMeshControl::SubStructure(), TPZBuildMultiphysicsMesh::TransferFromMeshes(), TPZBuildMultiphysicsMesh::TransferFromMultiPhysics(), TPZSubCompMesh::TransferMultiphysicsElementSolution(), TransferMultiphysicsSolution(), TPZMHMeshControl::TransferToMultiphysics(), UpdatePreviousState(), TPZHybridizeHDiv::VerifySolutionConsistency(), and TPZMHMeshControl::~TPZMHMeshControl().
|
inline |
|
inline |
Access method for the element solution vectors.
Definition at line 225 of file pzcmesh.h.
References AllocateNewConnect(), BuildCornerConnectList(), CleanUpUnconnectedNodes(), ComputeElGraph(), ComputeNodElCon(), ConnectSolution(), ExpandSolution(), fElementSolution, FindMaterial(), InitializeBlock(), InsertMaterialObject(), LoadReferences(), Print(), and SetElementSolution().
Referenced by TPZCompMeshTools::ComputeDifferenceNorm(), ComputeError(), TPZMultiphysicsCompEl< TGeometry >::EvaluateError(), TPZAnalysisError::EvaluateError(), TPZInterpolationSpace::EvaluateError(), TPZSubCompMesh::EvaluateError(), TPZCompEl::Solution(), and TPZPostProcAnalysis::TransferSolution().
|
inline |
Returns a reference to the element pointers vector.
Definition at line 183 of file pzcmesh.h.
References fElementVec.
Referenced by TPZAgglomerateElement::AccumulateVertices(), TPZMHMeshControl::AddBoundaryInterfaceElements(), TPZBuildMultiphysicsMesh::AddConnects(), TPZMultiphysicsCompMesh::AddConnects(), TPZElementGroup::AddElement(), TPZMHMeshControl::AddElementBoundaries(), TPZBuildMultiphysicsMesh::AddElements(), TPZMultiphysicsCompMesh::AddElements(), TPZAgglomerateElement::AddSubElementIndex(), TPZFrontStructMatrix< front >::AdjustSequenceNumbering(), TPZBuildMultiphysicsMesh::AppendConnects(), TPZFrontStructMatrix< front >::Assemble(), TPZFrontStructMatrix< front >::AssembleNew(), TPZHybridizeHDiv::AssociateElements(), TPZPostProcAnalysis::AutoBuildDisc(), TPZSubCompMesh::BuildCornerConnectList(), BuildCornerConnectList(), BuildElementGroups(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZCreateApproximationSpace::BuildMesh(), BuildTransferMatrixDesc(), Hdiv2dPaper201504::ChangeExternalOrderConnects(), hdivCurvedJCompAppMath::ChangeExternalOrderConnects(), Hdiv3dPaper201504::ChangeExternalOrderConnects(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), Coarsen(), TPZEulerAnalysis::CompareRhs(), TPZAgglomerateElement::ComputeNeighbours(), TPZCreateApproximationSpace::CondenseLocalEquations(), ConvertDiscontinuous2Continuous(), TPZDohrStructMatrix::CorrectNeighbourDomainIndex(), TPZAgglomerateElement::CreateAgglomerateMesh(), TCedricTest::CreateCondensedElements(), TPZMHMixedMeshControl::CreateHDivMHMMesh(), TPZCreateApproximationSpace::CreateInterfaceElements(), TPZCreateApproximationSpace::CreateInterfaces(), TPZMHMeshControl::CreateInternalElements(), TPZAnalysis::CreateListOfCompElsToComputeError(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CreateSkeleton(), DeltaX(), TPZCompElDisc::Divide(), TPZInterpolatedElement::Divide(), TPZParFrontStructMatrix< front >::ElementAssemble(), Hdiv2dPaper201504::ErrorH1(), hdivCurvedJCompAppMath::ErrorH1(), Hdiv3dPaper201504::ErrorH1(), Hdiv2dPaper201504::ErrorPrimalDual(), hdivCurvedJCompAppMath::ErrorPrimalDual(), Hdiv3dPaper201504::ErrorPrimalDual(), TPZAnalysisError::EvaluateError(), TPZCompElDisc::EvaluateSquareResidual2D(), TPZStructMatrixST::ExecuteAssemble(), TPZGraphMesh::FindFirstInterpolatedElement(), TPZInterfaceElement::FreeInterface(), TPZParFrontStructMatrix< front >::GlobalAssemble(), TPZSubCompMesh::HasMaterial(), TPZAnalysisError::HPAdapt(), TPZCreateApproximationSpace::Hybridize(), TPZDohrStructMatrix::IdentifyCornerNodes(), TPZElastoPlasticAnalysis::IdentifyEquationsToZero(), TCedricTest::InterpolationError(), LesserEdgeOfMesh(), TPZAgglomerateElement::ListOfGroupings(), TPZInterfaceElement::main(), TPZCreateApproximationSpace::MakeRaviartThomas(), MaximumRadiusOfMesh(), TPZFlowCompMesh::MaxVelocityOfMesh(), TPZMGAnalysis::MeshError(), TPZSubCompMesh::NeedsComputing(), TPZStructMatrixOR::ThreadData::NextElement(), TPZPairStructMatrix::ThreadData::NextElement(), TPZStructMatrixCS::ThreadData::NextElement(), NSubMesh(), TPZAnalysis::OptimizeBandwidth(), TPZFrontStructMatrix< front >::OrderElement(), TPZStructMatrixGCTP::OrderElement(), TPZStructMatrixGC::OrderElement(), TPZStructMatrixOT::OrderElement(), TPZAnalysis::PostProcess(), TPZAnalysis::PostProcessErrorSerial(), TPZAnalysis::PostProcessTable(), TPZAgglomerateElement::Print(), TPZMHMeshControl::PrintBoundaryInfo(), TPZVTKGeoMesh::PrintCMeshVTK(), TPZVTKGeoMesh::PrintPOrderPoints(), TPZMHMeshControl::PrintSubdomain(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), ProjectSolution(), TPZCompMeshReferred::ReferredEl(), RemakeAllInterfaceElements(), RemoveEl(), ResequenceByGeometry(), ResetMesh(), TPZNonLinMultGridAnalysis::ResetReference(), TPZGeoMesh::RestoreReference(), SaddlePermute(), SaddlePermute2(), TPZStructMatrixGCTP::Serial_Assemble(), TPZStructMatrixGC::Serial_Assemble(), TPZStructMatrixCS::Serial_Assemble(), TPZStructMatrixOT::Serial_Assemble(), TPZStructMatrixOR::Serial_Assemble(), TPZPairStructMatrix::SerialAssemble(), TPZGraphMesh::SetCompMesh(), TPZStructMatrixBase::SetMaterialIds(), TPZNonLinMultGridAnalysis::SetReference(), TPZCompElDisc::SetTensorialShape(), TPZCompElDisc::SetTotalOrderShape(), TPZAgglomerateElement::SubElement(), SubMesh(), TPZDohrStructMatrix::SubStructure(), SubStructure(), TPZGenSubStruct::SubStructure(), TPZPairStructMatrix::TBBAssemble(), TPZStructMatrixOR::ThreadData::ThreadWork(), TPZPairStructMatrix::ThreadData::ThreadWork(), TPZStructMatrixGC::ThreadData::ThreadWork(), TPZStructMatrixCS::ThreadData::ThreadWork(), TPZStructMatrixOT::ThreadData::ThreadWork(), TPZStructMatrixGC::ThreadData::ThreadWorkResidual(), TPZStructMatrixOT::ThreadData::ThreadWorkResidual(), TPZCompEl::TPZCompEl(), TPZCompElHDivBound2< TSHAPE >::TPZCompElHDivBound2(), TPZCondensedCompEl::TPZCondensedCompEl(), TPZGraphMesh::TPZGraphMesh(), TPZInterfaceElement::TPZInterfaceElement(), TPZMultiphysicsInterfaceElement::TPZMultiphysicsInterfaceElement(), TPZPairStructMatrix::TPZPairStructMatrix(), TPZSubCompMesh::TransferElementFrom(), TPZSubCompMesh::TransferElementTo(), TPZMHMeshControl::TransferToMultiphysics(), TPZCompMeshTools::UnCondensedElements(), TPZCreateApproximationSpace::UndoCondenseLocalEquations(), TPZCreateApproximationSpace::UndoMakeRaviartThomas(), TPZCompMeshTools::UnGroupElements(), TPZMGAnalysis::UniformlyRefineMesh(), TPZNonLinMultGridAnalysis::UniformlyRefineMesh(), TPZBuildMultiphysicsMesh::UniformRefineCompEl(), TPZBuildMultiphysicsMesh::UniformRefineCompMesh(), TPZCondensedCompEl::Unwrap(), TPZElementGroup::Unwrap(), TCedricTest::UnwrapElements(), TPZStructMatrixTBBFlow::Write(), TPZCompEl::~TPZCompEl(), and TPZCondensedCompEl::~TPZCondensedCompEl().
|
inline |
Returns a reference to the element pointers vector.
Definition at line 195 of file pzcmesh.h.
References fElementVec.
void TPZCompMesh::EvaluateError | ( | std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)> | fp, |
bool | store_error, | ||
TPZVec< REAL > & | errorSum | ||
) |
Evaluates the error given the two vectors of the analised parameters.
fp | pointer for the function with following arguments: |
errorSum | - return the L1 error |
Definition at line 1395 of file pzcmesh.cpp.
References TPZCompEl::EvaluateError(), fElementVec, TPZVec< T >::Fill(), TPZMaterial::Id(), TPZCompEl::Material(), TPZVec< T >::NElements(), TPZChunkVector< T, EXP >::NElements(), TPZVec< T >::Resize(), and sqrt.
Referenced by TPZAnalysisError::EvaluateError(), TPZAnalysisError::hp_Adaptive_Mesh_Design(), and SetAllCreateFunctionsContinuousWithMem().
|
virtual |
Adapt the solution vector to new block dimensions.
Definition at line 396 of file pzcmesh.cpp.
References TPZMatrix< TVar >::Cols(), TPZBlock< TVar >::Dim(), fBlock, fSolution, fSolutionBlock, TPZBlock< TVar >::NBlocks(), TPZBlock< TVar >::Position(), TPZFMatrix< TVar >::Redim(), TPZBlock< TVar >::Resequence(), and TPZBlock< TVar >::Size().
Referenced by TPZBuildMultiphysicsMesh::AddConnects(), TPZMultiphysicsCompMesh::AddConnects(), TPZCompMeshTools::AdjustFluxPolynomialOrders(), TPZBuildMultiphysicsMesh::AppendConnects(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedHybridMeshControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), Hdiv2dPaper201504::ChangeExternalOrderConnects(), hdivCurvedJCompAppMath::ChangeExternalOrderConnects(), Hdiv3dPaper201504::ChangeExternalOrderConnects(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), CreateCompMesh(), TPZMHMixedMeshControl::CreateHDivMHMMesh(), TPZMHMixedHybridMeshControl::CreateInternalFluxElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedHybridMeshControl::CreatePressureInterfaces(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedHybridMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CriaMalhaTemporaria(), Discontinuous2Continuous(), TPZCompElDisc::Divide(), TPZInterpolatedElement::Divide(), TPZCompMeshReferred::DivideReferredEl(), ElementSolution(), TPZFlowCompMesh::ExpandSolution2(), TPZCreateApproximationSpace::Hybridize(), TPZMHMixedMeshControl::HybridizeSkeleton(), TPZMHMixedHybridMeshControl::HybridizeSkeleton(), TPZMHMeshControl::HybridizeSkeleton(), InitializeBlock(), TPZSubCompMesh::MakeAllInternal(), TPZCreateApproximationSpace::MakeRaviartThomas(), Permute(), Hdiv2dPaper201504::PrintErrors(), hdivCurvedJCompAppMath::PrintErrors(), Hdiv3dPaper201504::PrintErrors(), Hdiv2dPaper201504::Run(), hdivCurvedJCompAppMath::Run(), Hdiv3dPaper201504::Run(), TPZSubCompMesh::SetNumberRigidBodyModes(), TPZPostProcAnalysis::SetPostProcessVariables(), TPZCompMeshTools::SetPressureOrders(), hdivCurvedJCompAppMath::SetupDisconnectedHdivboud(), TPZGenSubStruct::SubStructure(), TPZMHMeshControl::TransferToMultiphysics(), and TPZCreateApproximationSpace::UndoMakeRaviartThomas().
|
inlinevirtual |
Get the father meshes stack.
Reimplemented in TPZSubCompMesh.
Definition at line 337 of file pzcmesh.h.
Referenced by TPZSubCompMesh::CommonMesh(), CommonMesh(), TPZSubCompMesh::IsAllowedElement(), and TPZSubCompMesh::SetNumberRigidBodyModes().
TPZMaterial * TPZCompMesh::FindMaterial | ( | int | id | ) |
Find the material with identity id.
id | material id to be found |
Definition at line 297 of file pzcmesh.cpp.
References fMaterialVec.
Referenced by TPZPoroElastoPlasticAnalysis::AcceptSolution(), TPZElastoPlasticAnalysis::AcceptSolution(), TPZPostProcAnalysis::AutoBuildDisc(), TPZCreateApproximationSpace::BuildMesh(), TPZSBFemVolume::ComputeKMatrices(), TPZInterpolationSpace::ComputeNormal(), TPZSBFemVolume::ComputeSolution(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZBuildSBFem::CreateVolumetricElements(), TPZBuildSBFem::CreateVolumetricElementsFromSkeleton(), TPZDXGraphMesh::DrawMesh(), TPZVTKGraphMesh::DrawSolution(), TPZMVGraphMesh::DrawSolution(), TPZV3DGraphMesh::DrawSolution(), TPZDXGraphMesh::DrawSolution(), ElementSolution(), TPZGenSubStruct::GenerateMesh(), TPZMHMixedHybridMeshControl::InsertPeriferalHdivMaterialObjects(), TPZMHMixedHybridMeshControl::InsertPeriferalMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalMaterialObjects(), TPZHybridizeHDiv::InsertPeriferalMaterialObjects(), TPZMHMeshControl::InsertPeriferalMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalPressureMaterialObjects(), main(), TPZElastoPlasticAnalysis::ManageIterativeProcess(), TPZCompEl::Material(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), PYBIND11_MODULE(), TPZPoroElastoPlasticAnalysis::SetContributionTime(), TPZPoroElastoPlasticAnalysis::SetDeltaT(), TPZPostProcAnalysis::SetPostProcessVariables(), TPZSBFemVolume::Shape(), TPZAnalysis::ShowShape(), TPZSBFemVolume::Solution(), TPZIntelGen< TSHAPE >::TPZIntelGen(), TPZMeshSolution::TPZMeshSolution(), TPZMultiphysicsInterfaceElement::TPZMultiphysicsInterfaceElement(), TPZSubCompMesh::TransferElementFrom(), and TPZNonLinMultGridAnalysis::TwoGridAlgorithm().
|
inline |
Definition at line 398 of file pzcmesh.h.
References fDefaultOrder.
Referenced by Coarsen(), TPZInterpolatedElement::CreateMidSideConnect(), TPZNonLinMultGridAnalysis::SetDeltaTime(), TPZBuildMultiphysicsMesh::ShowShape(), TPZAnalysis::ShowShape(), TPZCompElDisc::TPZCompElDisc(), TPZCompElHDiv< TSHAPE >::TPZCompElHDiv(), TPZCompElHDivBound2< TSHAPE >::TPZCompElHDivBound2(), TPZCompElHDivPressure< TSHAPE >::TPZCompElHDivPressure(), and TPZInterpolationSpace::TPZInterpolationSpace().
void TPZCompMesh::GetElementPatch | ( | TPZVec< int64_t > | nodtoelgraph, |
TPZVec< int64_t > | nodtoelgraphindex, | ||
TPZStack< int64_t > & | elgraph, | ||
TPZVec< int64_t > & | elgraphindex, | ||
int64_t | elind, | ||
TPZStack< int64_t > & | patch | ||
) |
Gives the element patch.
Definition at line 1599 of file pzcmesh.cpp.
References TPZStack< T, NumExtAlloc >::Push(), and TPZManVector< T, NumExtAlloc >::Resize().
|
virtual |
Get an external connection from the supermesh - Supermesh is one mesh who contains the analised submesh.
superind | index of the element to be trasnfered |
super | pointer to the destination mesh |
Reimplemented in TPZSubCompMesh.
Definition at line 1517 of file pzcmesh.cpp.
Referenced by TPZInterfaceElement::ConnectIndex(), TPZInterfaceElement::GetConnects(), TPZSubCompMesh::GetFromSuperMesh(), and TransferElement().
|
inline |
Definition at line 449 of file pzcmesh.h.
References AutoBuild(), and fNmeshes.
Referenced by TPZBuildMultiphysicsMesh::AppendConnects().
void TPZCompMesh::GetNodeToElGraph | ( | TPZVec< int64_t > & | nodtoelgraph, |
TPZVec< int64_t > & | nodtoelgraphinde, | ||
TPZStack< int64_t > & | elgraph, | ||
TPZVec< int64_t > & | elgraphindexx | ||
) |
Gives the conects graphs.
Definition at line 1576 of file pzcmesh.cpp.
References ComputeElGraph(), TPZVec< T >::NElements(), NIndependentConnects(), and TPZRenumbering::SetElementGraph().
void TPZCompMesh::GetRefPatches | ( | std::set< TPZGeoEl *> & | grpatch | ) |
Gives all patches of the mesh.
grpatch | - stack where will be inserted the geometric elements |
Definition at line 1558 of file pzcmesh.cpp.
References fElementVec, and NElements().
void TPZCompMesh::InitializeBlock | ( | ) |
Resequence the block object, remove unconnected connect objects and reset the dimension of the solution vector.
Definition at line 391 of file pzcmesh.cpp.
References CleanUpUnconnectedNodes(), and ExpandSolution().
Referenced by AdjustBoundaryElements(), AutoBuildContDisc(), TPZPostProcAnalysis::AutoBuildDisc(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZCreateApproximationSpace::BuildMesh(), TPZDohrStructMatrix::Create(), TPZHybridizeHDiv::CreateInterfaceElements(), ElementSolution(), TPZMHMixedHybridMeshControl::GroupandCondenseElements(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZSymetricSpStructMatrix::main(), TPZSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), TPZSubCompMesh::main(), and TPZMHMeshControl::SubStructure().
int TPZCompMesh::InsertMaterialObject | ( | TPZMaterial * | mat | ) |
Insert a material object in the datastructure.
mat | pointer to the material |
Insert a material object in the datastructure
Definition at line 287 of file pzcmesh.cpp.
References DebugStop, fMaterialVec, and TPZMaterial::Id().
Referenced by TPZMulticamadaOrthotropic::AddPlacaOrtho(), Hdiv2dPaper201504::CMeshFlux(), hdivCurvedJCompAppMath::CMeshFlux(), Hdiv3dPaper201504::CMeshFlux(), Hdiv2dPaper201504::CMeshH1(), hdivCurvedJCompAppMath::CMeshH1(), Hdiv3dPaper201504::CMeshH1(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), CreatePressureMesh(), TPZMHMeshControl::CriaMalhaTemporaria(), ElementSolution(), TPZCompElDisc::EvaluateSquareResidual2D(), TCedricTest::GenerateCompMesh(), TPZGenSubStruct::GenerateMesh(), insert_elasticity(), InsertElasticity(), Input::InsertElasticityCubo(), InsertElasticityCubo(), TPZMHMixedHybridMeshControl::InsertPeriferalHdivMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalHdivMaterialObjects(), TPZMHMixedHybridMeshControl::InsertPeriferalMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalMaterialObjects(), TPZHybridizeHDiv::InsertPeriferalMaterialObjects(), TPZMHMeshControl::InsertPeriferalMaterialObjects(), TPZMHMixedHybridMeshControl::InsertPeriferalPressureMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalPressureMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalRotationMaterialObjects(), InsertViscoElasticity(), InsertViscoElasticityCubo(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZSymetricSpStructMatrix::main(), TPZSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), TPZSubCompMesh::main(), MalhaComp(), PYBIND11_MODULE(), TPZReadMeshHR::ReadBCs(), TPZReadMeshHR::ReadMaterials(), and TPZPostProcAnalysis::SetPostProcessVariables().
TPZVec< STATE > TPZCompMesh::Integrate | ( | const std::string & | varname, |
const std::set< int > & | matids | ||
) |
Integrate the variable name over the mesh.
Integrate the postprocessed variable name over the elements included in the set matids.
Definition at line 2162 of file pzcmesh.cpp.
References DebugStop, Element(), TPZCompEl::IntegrateSolution(), TPZBndCond::Material(), TPZGeoEl::MaterialId(), MaterialVec(), NElements(), TPZMaterial::NSolutionVariables(), TPZCompEl::Reference(), TPZVec< T >::size(), and TPZMaterial::VariableIndex().
Referenced by TPZSubCompMesh::IntegrateSolution(), and SetAllCreateFunctionsContinuousWithMem().
REAL TPZCompMesh::LesserEdgeOfMesh | ( | ) |
Definition at line 1849 of file pzcmesh.cpp.
References dist(), EAgglomerate, EInterface, ElementVec(), TPZMaterial::Id(), TPZCompEl::LesserEdgeOfEl(), TPZCompEl::Material(), TPZChunkVector< T, EXP >::NElements(), and TPZCompEl::Type().
Referenced by TPZFlowCompMesh::ComputeTimeStep(), SetAllCreateFunctionsContinuousWithMem(), and TPZNonLinMultGridAnalysis::SetDeltaTime().
void TPZCompMesh::LoadReferences | ( | ) |
Map this grid in the geometric grid.
Definition at line 482 of file pzcmesh.cpp.
References fElementVec, TPZCompEl::LoadElementReference(), NElements(), Reference(), and TPZGeoMesh::SetReference().
Referenced by TPZMHMeshControl::AddBoundaryInterfaceElements(), TPZMHMixedHybridMeshControl::AdjustBoundaryConditionsOfFractures(), TPZPostProcAnalysis::AutoBuildDisc(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), BuildElementGroups(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), BuildTransferMatrix(), BuildTransferMatrixDesc(), CleanUp(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), TPZHybridizeHDiv::CreateInterfaceElements(), TPZMHMeshControl::CreateInterfaceElements(), TPZMHMixedHybridMeshControl::CreateInternalFluxElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZHybridizeHDiv::CreateMultiphysicsMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZBuildSBFem::CreateVolumetricElements(), TPZBuildSBFem::CreateVolumetricElementsFromSkeleton(), TPZCompMeshReferred::DivideReferredEl(), ElementSolution(), TPZCompElDisc::EvaluateSquareResidual2D(), TPZCompMeshTools::GroupElements(), TPZMHMixedHybridMeshControl::GroupElements(), TPZMHMixedMeshChannelControl::HideTheElements(), TPZMHMixedMeshControl::HideTheElements(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZMHMixedMeshControl::HybridizeSkeleton(), TPZMHMeshControl::HybridizeSkeleton(), TPZMHMixedMeshControl::InsertPeriferalHdivMaterialObjects(), TPZSubCompMesh::LoadElementReference(), TPZCompMeshReferred::LoadReferred(), TPZMGAnalysis::MeshError(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), TPZVTKGeoMesh::PrintCMeshVTK(), ProjectSolution(), TPZCompMeshTools::SetPressureOrders(), TPZGenSubStruct::SubStructure(), TPZMHMeshControl::SubStructure(), TPZMHMeshControl::TransferToMultiphysics(), TPZHybridizeHDiv::VerifySolutionConsistency(), and TPZSubCompMesh::~TPZSubCompMesh().
void TPZCompMesh::LoadSolution | ( | const TPZFMatrix< STATE > & | sol | ) |
Given the solution of the global system of equations, computes and stores the solution for the restricted nodes.
sol | given solution matrix |
Definition at line 441 of file pzcmesh.cpp.
References TPZMatrix< TVar >::Cols(), fElementVec, fSolution, TPZFMatrix< TVar >::GetVal(), TPZCompEl::LoadSolution(), NElements(), TPZFMatrix< TVar >::Resize(), TPZMatrix< TVar >::Rows(), and val().
Referenced by TPZAnalysis::AnimateRun(), TPZMGAnalysis::AppendMesh(), TPZPlasticDiagnostic::CheckGlobal(), TPZEulerAnalysis::LineSearch(), TPZAnalysis::LoadSolution(), TPZSubCompMesh::LoadSolution(), main(), TPZAnalysis::PostProcess(), TPZAnalysis::PostProcessErrorParallel(), TPZAnalysis::PostProcessErrorSerial(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), TPZPoroElastoPlasticAnalysis::SetAdvancedState(), TPZEulerAnalysis::SetAdvancedState(), SetAllCreateFunctionsContinuousWithMem(), TPZPoroElastoPlasticAnalysis::SetLastState(), TPZEulerAnalysis::SetLastState(), TPZMGAnalysis::Solve(), TPZAnalysis::Solve(), and TPZEulerAnalysis::UpdateSolAndRhs().
|
inlinevirtual |
Make all mesh connections internal mesh connections. Connects to an internal connection.
Reimplemented in TPZSubCompMesh.
|
inlinevirtual |
Makes a specified connection a internal mesh connection.
local | connection local number to be processed |
Reimplemented in TPZSubCompMesh.
|
inline |
Returns a reference to the material pointers vector.
Definition at line 203 of file pzcmesh.h.
References fMaterialVec.
Referenced by TPZMGAnalysis::AppendMesh(), TPZCreateApproximationSpace::BuildMesh(), TPZMHMixedMeshControl::BuildMultiPhysicsMesh(), TPZSubCompMesh::CalcStiff(), TPZAnalysis::ComputeNumberofLoadCases(), TPZPoroElastoPlasticAnalysis::FindPorousMaterials(), TPZAnalysis::HighestDimension(), TPZAnalysis::IdentifyPostProcessingMatIds(), TPZSubCompMesh::InitializeEF(), TPZMHMixedHybridMeshControl::InsertPeriferalHdivMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalPressureMaterialObjects(), TPZMHMixedMeshControl::InsertPeriferalRotationMaterialObjects(), TPZAnalysis::Integrate(), Integrate(), main(), TPZMGAnalysis::MeshError(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), PYBIND11_MODULE(), TPZSubCompMesh::Read(), TPZReadMeshHR::ReadBCs(), TPZTransientAnalysis< TRANSIENTCLASS >::SetAllMaterialsDeltaT(), TPZTransientAnalysis< TRANSIENTCLASS >::SetCurrentState(), TPZTransientAnalysis< TRANSIENTCLASS >::SetExplicit(), TPZTransientAnalysis< TRANSIENTCLASS >::SetFluxOnly(), TPZTransientAnalysis< TRANSIENTCLASS >::SetImplicit(), TPZTransientAnalysis< TRANSIENTCLASS >::SetLastState(), TPZTransientAnalysis< TRANSIENTCLASS >::SetMassMatrix(), TPZElastoPlasticAnalysis::SetUpdateMem(), TPZAnalysis::ShowShape(), TPZCheckRestraint::TPZCheckRestraint(), TPZSubCompMesh::TransferElementFrom(), TPZMHMeshControl::TransferToMultiphysics(), TPZNonLinMultGridAnalysis::UniformlyRefineMesh(), and TPZSubCompMesh::~TPZSubCompMesh().
|
inline |
Returns a reference to the material pointers vector.
Definition at line 206 of file pzcmesh.h.
References fMaterialVec.
REAL TPZCompMesh::MaximumRadiusOfMesh | ( | ) |
Definition at line 1833 of file pzcmesh.cpp.
References dist(), EInterface, ElementVec(), TPZMaterial::Id(), TPZCompEl::Material(), TPZCompEl::MaximumRadiusOfEl(), TPZChunkVector< T, EXP >::NElements(), and TPZCompEl::Type().
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
protected |
Modify the permute vector swapping the lagrangeq with maxeq and shifting the intermediate equations.
Definition at line 2733 of file pzcmesh.cpp.
References DebugStop, input, and TPZVec< T >::size().
Referenced by SaddlePermute2().
|
inline |
|
inline |
Number of connects allocated including free nodes.
Definition at line 163 of file pzcmesh.h.
References TPZChunkVector< T, EXP >::NElements(), and NIndependentConnects().
Referenced by TPZBuildMultiphysicsMesh::AppendConnects(), TPZHybridizeHDiv::AssociateElements(), TPZPostProcAnalysis::AutoBuildDisc(), CleanUpUnconnectedNodes(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZRenumbering::ColorElements(), TPZBuildMultiphysicsMesh::ComputeAtomicIndexes(), ComputeNodElCon(), ConnectSolution(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZStructMatrixGC::ElementColoring(), TPZStructMatrixOT::ElementColoring(), TPZMHMeshControl::HybridizeSkeleton(), NEquations(), NIndependentConnects(), TPZStructMatrixGCTP::OrderElement(), TPZStructMatrixGC::OrderElement(), TPZStructMatrixOT::OrderElement(), Permute(), Print(), ResequenceByGeometry(), TPZHybridizeHDiv::SplitConnects(), TPZMHMeshControl::TransferToMultiphysics(), TPZBuildMultiphysicsMesh::UniformRefineCompEl(), TPZBuildMultiphysicsMesh::UniformRefineCompMesh(), and TPZMHMixedMeshControl::~TPZMHMixedMeshControl().
|
inline |
Number of computational elements allocated.
Definition at line 169 of file pzcmesh.h.
References TPZChunkVector< T, EXP >::NElements().
Referenced by TPZMHMeshControl::AddBoundaryInterfaceElements(), TPZBuildMultiphysicsMesh::AddConnects(), TPZMultiphysicsCompMesh::AddConnects(), TPZBuildMultiphysicsMesh::AddElements(), TPZMultiphysicsCompMesh::AddElements(), TPZCompMeshTools::AddHDivPyramidRestraints(), TPZCompMeshTools::AdjustFluxPolynomialOrders(), TPZBuildMultiphysicsMesh::AppendConnects(), TPZSubMeshAnalysis::Assemble(), TPZParFrontStructMatrix< front >::Assemble(), TPZFrontStructMatrix< front >::Assemble(), AssembleError(), TPZFrontStructMatrix< front >::AssembleNew(), TPZHybridizeHDiv::AssociateElements(), BandWidth(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedMeshControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), TPZSubCompMesh::BuildCornerConnectList(), BuildCornerConnectList(), BuildElementGroups(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZMHMixedMeshControl::BuildMultiPhysicsMesh(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory(), BuildTransferMatrix(), TPZMHMixedHybridMeshControl::CheckMeshConsistency(), TPZMultiphysicsCompMesh::CleanElementsConnects(), CleanUp(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZRenumbering::ColorElements(), TPZEulerAnalysis::CompareRhs(), TPZBuildMultiphysicsMesh::ComputeAtomicIndexes(), TPZCompMeshTools::ComputeDifferenceNorm(), ComputeElGraph(), TPZAgglomerateElement::ComputeNeighbours(), ComputeNodElCon(), TPZCreateApproximationSpace::CondenseLocalEquations(), ConvertDiscontinuous2Continuous(), TPZDohrStructMatrix::CorrectNeighbourDomainIndex(), TPZAgglomerateElement::CreateAgglomerateMesh(), CreateCompMesh(), TCedricTest::CreateCondensedElements(), TPZCompMeshTools::CreatedCondensedElements(), TPZBuildSBFem::CreateElementGroups(), TPZMHMixedHybridMeshControl::CreateHDivWrappers(), TPZHybridizeHDiv::CreateInterfaceElements(), TPZMHMixedHybridMeshControl::CreateInternalAxialFluxes(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedHybridMeshControl::CreateLowerDimensionPressureElements(), TPZMHMixedHybridMeshControl::CreateMultiPhysicsInterfaceElements(), TPZMHMixedMeshControl::CreateMultiPhysicsInterfaceElements(), TPZMHMixedHybridMeshControl::CreatePressureInterfaces(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedHybridMeshControl::CreateSkeletonAxialFluxes(), TPZMHMixedMeshControl::DeletePressureElements(), TPZStructMatrixGC::ElementColoring(), TPZStructMatrixOT::ElementColoring(), Hdiv2dPaper201504::ErrorH1(), hdivCurvedJCompAppMath::ErrorH1(), Hdiv3dPaper201504::ErrorH1(), Hdiv2dPaper201504::ErrorPrimalDual(), hdivCurvedJCompAppMath::ErrorPrimalDual(), Hdiv3dPaper201504::ErrorPrimalDual(), TPZCompElDisc::EvaluateSquareResidual2D(), TPZStructMatrixST::ExecuteAssemble(), TPZCompMeshTools::ExpandHDivPyramidRestraints(), TPZGraphMesh::FindFirstInterpolatedElement(), TPZInterfaceElement::FreeInterface(), GetRefPatches(), TPZHybridizeHDiv::GroupandCondenseElements(), TPZCompMeshTools::GroupElements(), TPZMHMixedHybridMeshControl::GroupElements(), TPZSubCompMesh::HasMaterial(), TPZMHMixedMeshChannelControl::HideTheElements(), TPZMHMixedMeshControl::HideTheElements(), TPZCreateApproximationSpace::Hybridize(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZDohrStructMatrix::IdentifyCornerNodes(), TPZElastoPlasticAnalysis::IdentifyEquationsToZero(), TPZAnalysis::Integrate(), Integrate(), TCedricTest::InterpolationError(), TPZAgglomerateElement::ListOfGroupings(), LoadReferences(), TPZCompMeshReferred::LoadReferred(), TPZCompMeshTools::LoadSolution(), LoadSolution(), main(), TPZInterfaceElement::main(), TPZCreateApproximationSpace::MakeRaviartThomas(), TPZMGAnalysis::MeshError(), TPZStructMatrixCS::MultiThread_Assemble(), TPZStructMatrixOT::MultiThread_Assemble(), NSubMesh(), TPZAnalysis::OptimizeBandwidth(), TPZStructMatrixGCTP::OrderElement(), TPZStructMatrixGC::OrderElement(), TPZStructMatrixOT::OrderElement(), Print(), TPZMHMeshControl::PrintBoundaryInfo(), TPZVTKGeoMesh::PrintCMeshVTK(), Hdiv2dPaper201504::PrintErrors(), hdivCurvedJCompAppMath::PrintErrors(), Hdiv3dPaper201504::PrintErrors(), TPZMHMixedHybridMeshControl::PrintFriendly(), TPZVTKGeoMesh::PrintPOrderPoints(), TPZMHMeshControl::PrintSubdomain(), TPZAnalysis::PrintVectorByElement(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), PYBIND11_MODULE(), ResetMesh(), TPZNonLinMultGridAnalysis::ResetReference(), TCedricTest::Run(), Hdiv2dPaper201504::Run(), hdivCurvedJCompAppMath::Run(), Hdiv3dPaper201504::Run(), SaddlePermute(), SaddlePermute2(), TPZStructMatrixGCTP::Serial_Assemble(), TPZStructMatrixGC::Serial_Assemble(), TPZStructMatrixCS::Serial_Assemble(), TPZStructMatrixOT::Serial_Assemble(), TPZStructMatrixOR::Serial_Assemble(), TPZPairStructMatrix::SerialAssemble(), SetElementSolution(), TPZCompMeshTools::SetPressureOrders(), TPZNonLinMultGridAnalysis::SetReference(), TPZCompElDisc::SetTensorialShape(), TPZCompElDisc::SetTotalOrderShape(), hdivCurvedJCompAppMath::SetupDisconnectedHdivboud(), Skyline(), TPZMHMixedHybridMeshControl::SplitFluxElementsAroundFractures(), SubMesh(), TPZDohrStructMatrix::SubStructure(), SubStructure(), TPZGenSubStruct::SubStructure(), TPZMHMeshControl::SubStructure(), TPZPairStructMatrix::TBBAssemble(), TPZStructMatrixOR::ThreadData::ThreadAssembly(), TPZPairStructMatrix::ThreadData::ThreadAssembly1(), TPZPairStructMatrix::ThreadData::ThreadAssembly2(), TPZStructMatrixOR::ThreadData::ThreadWork(), TPZPairStructMatrix::ThreadData::ThreadWork(), TPZStructMatrixGC::ThreadData::ThreadWork(), TPZStructMatrixCS::ThreadData::ThreadWork(), TPZStructMatrixGC::ThreadData::ThreadWorkResidual(), TPZPairStructMatrix::TPZPairStructMatrix(), TPZSubCompMesh::TransferElementTo(), TPZBuildMultiphysicsMesh::TransferFromMeshes(), TPZBuildMultiphysicsMesh::TransferFromMultiPhysics(), TPZSubCompMesh::TransferMultiphysicsElementSolution(), TransferMultiphysicsSolution(), TPZPostProcAnalysis::TransferSolution(), TPZMHMeshControl::TransferToMultiphysics(), TPZCompMeshTools::UnCondensedElements(), TPZCreateApproximationSpace::UndoCondenseLocalEquations(), TPZCreateApproximationSpace::UndoMakeRaviartThomas(), TPZCompMeshTools::UnGroupElements(), TPZBuildMultiphysicsMesh::UniformRefineCompEl(), TPZBuildMultiphysicsMesh::UniformRefineCompMesh(), UpdatePreviousState(), TPZHybridizeHDiv::VerifySolutionConsistency(), TPZMHMeshControl::~TPZMHMeshControl(), and TPZSubCompMesh::~TPZSubCompMesh().
int64_t TPZCompMesh::NEquations | ( | ) |
This computes the number of equations associated with non-restrained nodes.
Definition at line 721 of file pzcmesh.cpp.
References DebugStop, fBlock, fConnectVec, TPZConnect::HasDependency(), TPZConnect::IsCondensed(), NConnects(), TPZConnect::NElConnected(), TPZConnect::NShape(), TPZConnect::NState(), TPZConnect::SequenceNumber(), and TPZBlock< TVar >::Size().
Referenced by TPZAnalysis::AnimateRun(), TPZMGAnalysis::AppendMesh(), TPZTransientAnalysis< TRANSIENTCLASS >::Assemble(), TPZSubMeshAnalysis::Assemble(), TPZDohrStructMatrix::Assemble(), TPZParFrontStructMatrix< front >::Assemble(), TPZFrontStructMatrix< front >::Assemble(), TPZAnalysis::Assemble(), TPZFrontStructMatrix< front >::AssembleNew(), TPZAnalysis::AssembleResidual(), TPZEulerAnalysis::BufferLastStateAssemble(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedHybridMeshControl::BuildComputationalMesh(), TPZMHMixedMeshControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), TPZAnalysis::BuildPreconditioner(), TPZSubCompMesh::CalcStiff(), TPZElastoPlasticAnalysis::CheckConv(), TPZPlasticDiagnostic::CheckGlobal(), ComputeFillIn(), TPZTransientAnalysis< TRANSIENTCLASS >::ComputeFluxOnly(), TPZTransientAnalysis< TRANSIENTCLASS >::ComputeLinearTangentMatrix(), TPZTransientAnalysis< TRANSIENTCLASS >::ComputeMassMatrix(), TPZNonLinearAnalysis::ComputeTangent(), TPZElastoPlasticAnalysis::ComputeTangent(), TPZMatRedStructMatrix< TStructMatrix, TSparseMatrix >::Create(), TPZDohrStructMatrix::Create(), TPZAnalysis::CreateListOfCompElsToComputeError(), TPZEulerAnalysis::EvaluateFluxEpsilon(), TPZGenSubStruct::GenerateMesh(), TPZElastoPlasticAnalysis::GetActiveEquations(), TPZAnalysisError::hp_Adaptive_Mesh_Design(), TPZGenSubStruct::IdentifyEqNumbers(), TPZGenSubStruct::InitializeDohr(), TPZGenSubStruct::InitializeDohrCondense(), TCedricTest::InterpolationError(), TPZNonLinearAnalysis::IterativeProcess(), TPZElastoPlasticAnalysis::IterativeProcess(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), main(), TPZAnalysis::PostProcessErrorSerial(), Hdiv2dPaper201504::PrintErrors(), hdivCurvedJCompAppMath::PrintErrors(), Hdiv3dPaper201504::PrintErrors(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), ProjectSolution(), PYBIND11_MODULE(), TPZNonLinearAnalysis::Residual(), TPZElastoPlasticAnalysis::Residual(), TCedricTest::Run(), TPZAnalysis::Run(), TPZTransientAnalysis< TRANSIENTCLASS >::RunTransient(), TPZSubMeshAnalysis::SetCompMesh(), TPZAnalysis::SetCompMesh(), SetDefaultOrder(), TPZStructMatrixBase::SetMesh(), Skyline(), TPZMGAnalysis::Solve(), TPZEulerAnalysis::Solve(), TPZAnalysis::Solve(), Hdiv2dPaper201504::SolveSyst(), hdivCurvedJCompAppMath::SolveSyst(), Hdiv3dPaper201504::SolveSyst(), TPZElastoPlasticAnalysis::TPZElastoPlasticAnalysis(), TPZPoroElastoPlasticAnalysis::TPZPoroElastoPlasticAnalysis(), and TPZSubMeshAnalysis::TPZSubMeshAnalysis().
int64_t TPZCompMesh::NIndependentConnects | ( | ) |
Number of independent connect objects.
Definition at line 1080 of file pzcmesh.cpp.
References fConnectVec, TPZConnect::HasDependency(), TPZConnect::IsCondensed(), NConnects(), and TPZConnect::SequenceNumber().
Referenced by TPZAnalysis::BuildPreconditioner(), BuildTransferMatrix(), BuildTransferMatrixDesc(), ComputeElGraph(), TPZDohrStructMatrix::Create(), GetNodeToElGraph(), TPZGenSubStruct::IdentifyCornerNodes(), TPZDohrStructMatrix::IdentifyCornerNodes(), NConnects(), TPZAnalysis::OptimizeBandwidth(), SaddlePermute(), SaddlePermute2(), and TPZDohrStructMatrix::SubStructure().
|
inline |
Number of materials.
Definition at line 172 of file pzcmesh.h.
Referenced by BuildTransferMatrix(), BuildTransferMatrixDesc(), CleanUp(), TPZPoroElastoPlasticAnalysis::FindPorousMaterials(), Print(), ProjectSolution(), and PYBIND11_MODULE().
TPZCompMesh & TPZCompMesh::operator= | ( | const TPZCompMesh & | copy | ) |
copy the content of the mesh
Definition at line 1683 of file pzcmesh.cpp.
References CleanUp(), TPZCompEl::Clone(), CopyMaterials(), fBlock, fConnectVec, fCreate, fDefaultOrder, fDimModel, fElementSolution, fElementVec, fName, fReference, fSolution, fSolutionBlock, TPZChunkVector< T, EXP >::NElements(), TPZGeoMesh::ResetReference(), TPZAdmChunkVector< T, EXP >::Resize(), and TPZBlock< TVar >::SetMatrix().
Referenced by TPZMultiphysicsCompMesh::operator=().
void TPZCompMesh::Permute | ( | TPZVec< int64_t > & | permute | ) |
Permute the sequence number of the connect objects It is a permute gather operation.
permute | vector of elements to permute |
ExpandSolution must be called before calling this
Definition at line 1321 of file pzcmesh.cpp.
References ExpandSolution(), fBlock, fConnectVec, fSolution, fSolutionBlock, TPZCompEl::Index(), LOGPZ_DEBUG, TPZBlock< TVar >::NBlocks(), NConnects(), TPZVec< T >::NElements(), TPZBlock< TVar >::Position(), TPZBlock< TVar >::Resequence(), TPZMatrix< TVar >::Rows(), TPZConnect::SequenceNumber(), TPZBlock< TVar >::Set(), TPZConnect::SetSequenceNumber(), TPZBlock< TVar >::Size(), and TPZVec< T >::size().
Referenced by TPZFrontStructMatrix< front >::AdjustSequenceNumbering(), AssembleMatrices(), CleanUpUnconnectedNodes(), TPZDohrStructMatrix::Create(), InitializeMatrices(), TPZAnalysis::OptimizeBandwidth(), TPZSubCompMesh::PermuteExternalConnects(), TPZSubCompMesh::PermuteInternalFirst(), ResequenceByGeometry(), SaddlePermute(), SaddlePermute2(), and SetAllCreateFunctionsContinuousWithMem().
|
virtual |
Prints mesh data.
out | indicates the device where the data will be printed |
Reimplemented in TPZSubCompMesh, and TPZCompMeshReferred.
Definition at line 236 of file pzcmesh.cpp.
References DebugStop, Dimension(), fConnectVec, fElementVec, fMaterialVec, fName, TPZGeoEl::Index(), NConnects(), NElements(), NMaterials(), TPZMaterial::Print(), TPZCompEl::Print(), and TPZCompEl::Reference().
Referenced by TPZCompMeshTools::AddHDivPyramidRestraints(), TPZParFrontStructMatrix< front >::Assemble(), TPZPostProcAnalysis::AutoBuildDisc(), TPZMHMixedHybridMeshControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), TPZMHMixedMeshControl::CreateHDivMHMMesh(), TPZMHMixedMeshControl::CreateHDivPressureMHMMesh(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZCompElDisc::Divide(), ElementSolution(), TPZGenSubStruct::GenerateMesh(), TPZParSkylineStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), main(), TPZCompMeshReferred::Print(), TPZMHMixedMeshControl::Print(), TPZSubCompMesh::Print(), TPZMHMeshControl::Print(), PYBIND11_MODULE(), TCedricTest::Run(), TPZMHMeshControl::TransferToMultiphysics(), TPZNonLinMultGridAnalysis::TwoGridAlgorithm(), and ~TPZCompMesh().
void TPZCompMesh::ProjectSolution | ( | TPZFMatrix< STATE > & | projectsol | ) |
Definition at line 1923 of file pzcmesh.cpp.
References DebugStop, TPZMaterial::Dimension(), TPZCompEl::Dimension(), EAgglomerate, ElementVec(), fMaterialVec, LoadReferences(), TPZChunkVector< T, EXP >::NElements(), NEquations(), NMaterials(), TPZAgglomerateElement::ProjectSolution(), PZError, TPZFMatrix< TVar >::Redim(), Reference(), TPZGeoMesh::ResetReference(), TPZCompEl::Type(), and TPZFMatrix< TVar >::Zero().
Referenced by SetDefaultOrder().
|
virtual |
Put an local connection in the supermesh - Supermesh is one mesh who contains the analised submesh.
local | local index of the element to be trasnfered |
super | pointer to the destination mesh |
Reimplemented in TPZSubCompMesh.
Definition at line 1512 of file pzcmesh.cpp.
Referenced by TPZInterfaceElement::ConnectIndex(), TPZInterfaceElement::GetConnects(), TPZSubCompMesh::PutinSuperMesh(), and TransferElement().
|
overridevirtual |
Read the element data from a stream.
Read the element data from a stream
Reimplemented from TPZSavable.
Reimplemented in TPZSubCompMesh, TPZFlowCompMesh, and TPZCompMeshReferred.
Definition at line 2006 of file pzcmesh.cpp.
References fBlock, fConnectVec, fCreate, fDefaultOrder, fDimModel, fElementSolution, fElementVec, fGMesh, fMaterialVec, fName, fNmeshes, fReference, fSolution, fSolutionBlock, TPZPersistenceManager::GetAutoPointer(), TPZPersistenceManager::GetInstance(), TPZCreateApproximationSpace::Read(), TPZAdmChunkVector< T, EXP >::Read(), TPZStream::Read(), TPZBlock< TVar >::Read(), TPZFMatrix< TVar >::Read(), and TPZStream::ReadPointers().
Referenced by main(), TPZCompMeshReferred::Read(), TPZFlowCompMesh::Read(), TPZSubCompMesh::Read(), and SetAllCreateFunctionsContinuousWithMem().
|
inline |
Returns a pointer to the geometrical mesh associated.
Definition at line 209 of file pzcmesh.h.
References fReference.
Referenced by TPZBuildMultiphysicsMesh::AddElements(), TPZMultiphysicsCompMesh::AddElements(), AutoBuild(), AutoBuildContDisc(), TPZPostProcAnalysis::AutoBuildDisc(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZBuildSBFem::BuildComputationalMeshFromSkeleton(), TPZBuildSBFem::BuildComputationMesh(), TPZBuildMultiphysicsMesh::BuildHybridMesh(), TPZCreateApproximationSpace::BuildMesh(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpace(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory(), BuildTransferMatrix(), BuildTransferMatrixDesc(), CleanUp(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), TPZMultiphysicsElement::CreateInterface(), TPZInterpolationSpace::CreateInterface(), TPZHybridizeHDiv::CreateInterfaceElements(), TPZMHMixedMeshControl::CreateMultiPhysicsInterfaceElements(), TPZHybridizeHDiv::CreateMultiphysicsMesh(), TPZBuildSBFem::CreateVolumetricElements(), TPZBuildSBFem::CreateVolumetricElementsFromSkeleton(), TPZCompElDisc::Divide(), TPZCompMeshReferred::DivideReferredEl(), TPZCompElDisc::EvaluateSquareResidual2D(), TPZMeshSolution::Execute(), TPZMultiphysicsCompEl< TGeometry >::GetReferenceIndexVec(), TPZCompMeshTools::GroupElements(), TPZMHMixedHybridMeshControl::GroupElements(), TPZMHMixedMeshChannelControl::HideTheElements(), TPZMHMixedMeshControl::HideTheElements(), TPZCreateApproximationSpace::Hybridize(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZDohrStructMatrix::IdentifyCornerNodes(), TCedricTest::LoadInterpolation(), LoadReferences(), TPZCompMeshReferred::LoadReferred(), main(), TPZAnalysisError::MathematicaPlot(), TPZMGAnalysis::MeshError(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), TPZAnalysis::PostProcessErrorSerial(), TPZInterpolatedElement::PRefine(), TPZVTKGeoMesh::PrintCMeshVTK(), TPZGradientReconstruction::ProjectionL2GradientReconstructed(), ProjectSolution(), TPZSubCompMesh::Read(), TPZCompEl::Reference(), TPZPostProcAnalysis::SetCompMesh(), TPZAnalysis::SetCompMesh(), TPZCompMeshTools::SetPressureOrders(), TPZCompElDisc::SizeOfElement(), TPZDohrStructMatrix::SubStructure(), SubStructure(), TPZGenSubStruct::SubStructure(), TPZMeshSolution::TPZMeshSolution(), TPZNonLinMultGridAnalysis::TwoGridAlgorithm(), TPZMGAnalysis::UniformlyRefineMesh(), TPZNonLinMultGridAnalysis::UniformlyRefineMesh(), TPZGenSubStruct::UniformRefine(), TPZHybridizeHDiv::VerifySolutionConsistency(), ~TPZCompMesh(), and TPZSubCompMesh::~TPZSubCompMesh().
void TPZCompMesh::RemakeAllInterfaceElements | ( | ) |
Deletes all interfaces and rebuild them all.
Definition at line 1281 of file pzcmesh.cpp.
References TPZInterpolationSpace::CreateInterfaces(), ElementVec(), TPZChunkVector< T, EXP >::NElements(), PZError, and TPZInterpolationSpace::RemoveInterfaces().
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
inlinevirtual |
Returns the rootmesh who have the specified connection.
Reimplemented in TPZSubCompMesh.
void TPZCompMesh::SaddlePermute | ( | ) |
Put the sequence number of the pressure connects after the seq number of the flux connects.
Definition at line 2328 of file pzcmesh.cpp.
References TPZCompEl::Connect(), ConnectVec(), DebugStop, ElementVec(), TPZConnect::HasDependency(), TPZConnect::IsCondensed(), TPZConnect::LagrangeMultiplier(), LOGPZ_DEBUG, Max(), Min(), TPZCompEl::NConnects(), NElements(), TPZChunkVector< T, EXP >::NElements(), NIndependentConnects(), Permute(), TPZConnect::Print(), TPZCompEl::Print(), TPZVec< T >::Resize(), TPZConnect::SequenceNumber(), TPZVec< T >::size(), and switchEq().
Referenced by TPZAnalysis::OptimizeBandwidth(), SetAllCreateFunctionsContinuousWithMem(), TPZSubCompMesh::SetAnalysisSkyline(), and TPZMHMeshControl::SubStructure().
void TPZCompMesh::SaddlePermute2 | ( | ) |
Definition at line 2588 of file pzcmesh.cpp.
References TPZCompEl::Connect(), ElementVec(), TPZConnect::LagrangeMultiplier(), LOGPZ_DEBUG, ModifyPermute(), TPZCompEl::NConnects(), NElements(), NIndependentConnects(), Permute(), TPZConnect::Print(), TPZVec< T >::Resize(), TPZConnect::SequenceNumber(), and TPZVec< T >::size().
Referenced by SetAllCreateFunctionsContinuousWithMem().
|
inline |
Definition at line 537 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctions().
Referenced by TPZMHMixedHybridMeshControl::CreateAxialFluxElement(), TPZCompElDisc::Divide(), and TPZInterpolatedElement::Divide().
|
inline |
Definition at line 508 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctionsContinuous().
Referenced by TPZBuildMultiphysicsMesh::BuildHybridMesh(), Hdiv2dPaper201504::CMeshH1(), hdivCurvedJCompAppMath::CMeshH1(), Hdiv3dPaper201504::CMeshH1(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZMHMeshControl::CreateInternalElements(), TPZMHMixedMeshControl::HybridizeSkeleton(), PYBIND11_MODULE(), TPZInterpolatedElement::SetCreateFunctions(), and TPZCompEl::SetCreateFunctions().
|
inline |
Definition at line 518 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctionsContinuousReferred().
Referenced by TPZReferredCompEl< TCOMPEL >::SetCreateFunctions().
|
inline |
Definition at line 552 of file pzcmesh.h.
References AdjustBoundaryElements(), Check(), ClassId(), Clone(), Coarsen(), ComputeFillIn(), ComputeMesh(), Consolidate(), ConvertDiscontinuous2Continuous(), CopyMaterials(), DeltaX(), Discontinuous2Continuous(), Divide(), EvaluateError(), Integrate(), LesserEdgeOfMesh(), LoadSolution(), MaximumRadiusOfMesh(), Permute(), Read(), RemakeAllInterfaceElements(), SaddlePermute(), SaddlePermute2(), TPZCreateApproximationSpace::SetAllCreateFunctionsContinuousWithMem(), TransferMultiphysicsSolution(), UpdatePreviousState(), val(), and Write().
Referenced by PYBIND11_MODULE(), and TPZCompElWithMem< TBASE >::SetCreateFunctions().
|
inline |
Definition at line 503 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctionsDiscontinuous().
Referenced by Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMeshControl::CriaMalhaTemporaria(), and TPZCompElDisc::SetCreateFunctions().
|
inline |
Definition at line 513 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctionsDiscontinuousReferred().
Referenced by TPZReferredCompEl< TCOMPEL >::SetCreateFunctions().
|
inline |
Definition at line 523 of file pzcmesh.h.
References Dimension(), and TPZCreateApproximationSpace::SetAllCreateFunctionsHDiv().
Referenced by Hdiv2dPaper201504::CMeshFlux(), hdivCurvedJCompAppMath::CMeshFlux(), Hdiv3dPaper201504::CMeshFlux(), TPZMHMixedHybridMeshControl::CreateAxialFluxElement(), TPZMHMixedMeshControl::HybridizeSkeleton(), TPZCompElHDiv< TSHAPE >::SetCreateFunctions(), and TPZCompElHDivBound2< TSHAPE >::SetCreateFunctions().
|
inline |
Definition at line 531 of file pzcmesh.h.
References Dimension(), and TPZCreateApproximationSpace::SetAllCreateFunctionsHDivPressure().
Referenced by MalhaComp(), TPZCompElHDivPressure< TSHAPE >::SetCreateFunctions(), and TPZCompElHDivPressureBound< TSHAPE >::SetCreateFunctions().
|
inline |
Definition at line 542 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctionsMultiphysicElem().
Referenced by TPZMHMixedMeshControl::BuildMultiPhysicsMesh(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpace(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), TPZMHMixedMeshControl::CreateHDivPressureMHMMesh(), TPZMultiphysicsCompEl< TGeometry >::SetCreateFunctions(), and TPZMHMeshControl::TransferToMultiphysics().
|
inline |
Definition at line 547 of file pzcmesh.h.
References TPZCreateApproximationSpace::SetAllCreateFunctionsMultiphysicElemWithMem().
Referenced by TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory().
|
inline |
Definition at line 403 of file pzcmesh.h.
References AssembleError(), BandWidth(), BuildTransferMatrix(), BuildTransferMatrixDesc(), NEquations(), ProjectSolution(), and Skyline().
Referenced by TPZBuildMultiphysicsMesh::AddWrap(), Hdiv2dPaper201504::CMeshFlux(), hdivCurvedJCompAppMath::CMeshFlux(), Hdiv3dPaper201504::CMeshFlux(), Hdiv2dPaper201504::CMeshH1(), hdivCurvedJCompAppMath::CMeshH1(), Hdiv3dPaper201504::CMeshH1(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZMHMixedHybridMeshControl::CreateAxialFluxElement(), Input::CreateCuboSkyMatrix(), TPZMHMixedHybridMeshControl::CreateHDivWrappers(), TPZMHMixedHybridMeshControl::CreateInternalFluxElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), TPZMHMixedHybridMeshControl::CreatePressureInterfaces(), CreatePressureMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CreateSkeleton(), TPZInterpolatedElement::Divide(), TPZMulticamadaOrthotropic::GenerateMesh(), TPZCreateApproximationSpace::Hybridize(), TPZMHMixedMeshControl::HybridizeSkeleton(), TPZMHMeshControl::HybridizeSkeleton(), TPZMHMixedMeshControl::InsertPeriferalHdivMaterialObjects(), TPZParSkylineStructMatrix::main(), TPBSpStructMatrix::main(), TPZSpStructMatrix::main(), TPZSymetricSpStructMatrix::main(), TPZFrontStructMatrix< front >::main(), TPZParFrontStructMatrix< front >::main(), main(), PYBIND11_MODULE(), TPZMHMeshControl::SetInternalPOrder(), and TPZHybridizeHDiv::SplitConnects().
|
inline |
Set de dimension of the domain of the problem.
Definition at line 139 of file pzcmesh.h.
References SetReference().
Referenced by TPZMHMeshControl::BuildComputationalMesh(), Hdiv2dPaper201504::CMeshFlux(), hdivCurvedJCompAppMath::CMeshFlux(), Hdiv3dPaper201504::CMeshFlux(), Hdiv2dPaper201504::CMeshH1(), hdivCurvedJCompAppMath::CMeshH1(), Hdiv3dPaper201504::CMeshH1(), Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), Hdiv3dPaper201504::CMeshMixed(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), TPZMHMixedHybridMeshControl::CreateAxialFluxElement(), Input::CreateCuboSkyMatrix(), TPZMHMixedMeshControl::CreateHDivPressureMHMMesh(), TPZMHMixedHybridMeshControl::CreateInternalFluxElements(), TPZMHMixedMeshControl::CreateInternalFluxElements(), TPZMHMeshControl::CreateLagrangeMultiplierMesh(), CreatePressureMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZMHMixedMeshControl::CreateSkeleton(), TPZMHMeshControl::CreateSkeleton(), TPZMHMeshControl::CriaMalhaTemporaria(), TCedricTest::GenerateCompMesh(), TPZHybridizeHDiv::HybridizeInternalSides(), TPZMHMeshControl::HybridizeSkeleton(), insert_elasticity(), InsertElasticity(), Input::InsertElasticityCubo(), InsertElasticityCubo(), TPZMHMixedMeshControl::InsertPeriferalHdivMaterialObjects(), InsertViscoElasticity(), InsertViscoElasticityCubo(), TPZAgglomerateElement::ListOfGroupings(), main(), MalhaComp(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), PYBIND11_MODULE(), TPZCompMesh(), TPZMHMeshControl::TPZMHMeshControl(), TPZMHMixedMeshControl::TPZMHMixedMeshControl(), TPZSubCompMesh::TPZSubCompMesh(), TPZMHMeshControl::TransferToMultiphysics(), and TPZNonLinMultGridAnalysis::TwoGridAlgorithm().
void TPZCompMesh::SetElementSolution | ( | int64_t | i, |
TPZVec< STATE > & | sol | ||
) |
Set a ith element solution, expanding the element-solution matrix if necessary.
Definition at line 1533 of file pzcmesh.cpp.
References TPZMatrix< TVar >::Cols(), fElementSolution, LOGPZ_DEBUG, NElements(), TPZVec< T >::NElements(), TPZFMatrix< TVar >::Resize(), TPZVec< T >::size(), and sqrt.
Referenced by ElementSolution().
void TPZCompMesh::SetName | ( | const std::string & | nm | ) |
Set the mesh name.
Definition at line 232 of file pzcmesh.cpp.
References fName.
Referenced by TPZMHMixedMeshControl::CreateHDivMHMMesh(), TPZMHMixedMeshControl::CreatePressureMHMMesh(), TPZMHMixedMeshControl::CreateRotationMesh(), TPZSubCompMesh::main(), TPZNonLinMultGridAnalysis::OneGridAlgorithm(), TPZCompMesh(), TPZNonLinMultGridAnalysis::TPZNonLinMultGridAnalysis(), and TPZNonLinMultGridAnalysis::TwoGridAlgorithm().
|
inline |
Definition at line 445 of file pzcmesh.h.
Referenced by TPZBuildMultiphysicsMesh::AddConnects(), TPZMultiphysicsCompMesh::BuildMultiphysicsSpace(), and TPZMultiphysicsCompMesh::BuildMultiphysicsSpaceWithMemory().
|
inline |
Sets the geometric reference mesh.
Definition at line 727 of file pzcmesh.h.
References TPZGeoMesh::Dimension(), fDimModel, fGMesh, and fReference.
Referenced by TPZMHMeshControl::operator=(), and SetDimModel().
|
inline |
Sets the geometric reference mesh.
Definition at line 733 of file pzcmesh.h.
References TPZGeoMesh::Dimension(), fDimModel, fGMesh, and fReference.
|
virtual |
This method computes the skyline of the system of equations.
skyline | vector where the skyline will be computed |
Definition at line 787 of file pzcmesh.cpp.
References TPZCompEl::BuildConnectList(), ConnectVec(), DebugStop, fBlock, fConnectVec, fElementVec, TPZConnect::HasDependency(), TPZConnect::IsCondensed(), NElements(), TPZVec< T >::NElements(), TPZChunkVector< T, EXP >::NElements(), NEquations(), TPZBlock< TVar >::Position(), TPZManVector< T, NumExtAlloc >::Resize(), TPZVec< T >::Resize(), TPZConnect::SequenceNumber(), TPZBlock< TVar >::Size(), and TPZVec< T >::size().
Referenced by TPZParSkylineStructMatrix::Create(), TPZSkylineStructMatrix::Create(), SetDefaultOrder(), and TPZSubCompMesh::SkylineInternal().
|
inline |
Access the solution vector.
Definition at line 219 of file pzcmesh.h.
References fSolution.
Referenced by TPZInterpolatedElement::AdjustPreferredSideOrder(), TPZMHMixedMeshChannelControl::BuildComputationalMesh(), TPZMHMixedHybridMeshControl::BuildComputationalMesh(), TPZMHMixedMeshControl::BuildComputationalMesh(), TPZMHMeshControl::BuildComputationalMesh(), TPZInterpolatedElement::CalcIntegral(), TPZAgglomerateElement::CalcStiff(), TPZPlasticDiagnostic::CheckGlobal(), TPZReducedSpace::ComputeSolution(), TPZInterpolatedElement::ComputeSolution(), TPZCompElDisc::ComputeSolution(), TPZCompElHDiv< TSHAPE >::ComputeSolutionHDiv(), TPZCompElHDivPressure< TSHAPE >::ComputeSolutionPressureHDiv(), ConnectSolution(), ConvertDiscontinuous2Continuous(), TPZMGAnalysis::ElementError(), TCedricTest::InterpolationError(), TPZElastoPlasticAnalysis::IterativeProcess(), TCedricTest::LoadInterpolation(), TPZSBFemElementGroup::LoadSolution(), TPZCompEl::LoadSolution(), TPZSubCompMesh::LoadSolution(), TPZMultiphysicsCompMesh::LoadSolutionFromMultiPhysics(), TPZAnalysisError::MathematicaPlot(), PYBIND11_MODULE(), TPZAnalysis::SetCompMesh(), TPZTransientAnalysis< TRANSIENTCLASS >::SetInitialSolution(), TPZTransientAnalysis< TRANSIENTCLASS >::SetInitialSolutionAsZero(), TPZCompElDisc::SolutionX(), TPZNonLinearAnalysis::TPZNonLinearAnalysis(), TPZBuildMultiphysicsMesh::TransferFromMeshes(), TPZBuildMultiphysicsMesh::TransferFromMultiPhysics(), and TPZMultiphysicsElement::TransferMultiphysicsElementSolution().
|
inline |
|
inlinevirtual |
Transfer one element form a submesh to another mesh.
mesh | mesh pointer who will receive the element |
elindex | element index to transfer |
Reimplemented in TPZSubCompMesh.
Definition at line 370 of file pzcmesh.h.
References CommonMesh(), GetFromSuperMesh(), and PutinSuperMesh().
|
inlinevirtual |
Transfer one element from a specified mesh to the current submesh.
mesh | pointer to the mesh whose the element from |
elindex | element index to transfer |
Reimplemented in TPZSubCompMesh.
Definition at line 356 of file pzcmesh.h.
Referenced by TPZSubCompMesh::TransferElementFrom().
|
inlinevirtual |
Transfer one element from a submesh to another mesh.
mesh | mesh pointer who will receive the element |
elindex | element index to transfer |
Reimplemented in TPZSubCompMesh.
Definition at line 363 of file pzcmesh.h.
Referenced by TPZSubCompMesh::TransferElement(), and TPZSubCompMesh::TransferElementTo().
void TPZCompMesh::TransferMultiphysicsSolution | ( | ) |
Transfer multiphysics mesh solution.
Definition at line 470 of file pzcmesh.cpp.
References Element(), NElements(), and TPZCompEl::TransferMultiphysicsElementSolution().
Referenced by SetAllCreateFunctionsContinuousWithMem(), TPZAnalysis::ShowShape(), and TPZAnalysis::Solve().
void TPZCompMesh::UpdatePreviousState | ( | STATE | mult = 1.0 | ) |
update the solution at the previous state with fSolution and set fSolution to the previous state
Definition at line 2826 of file pzcmesh.cpp.
References TPZMatrix< TVar >::Cols(), Element(), fSolN, fSolution, NElements(), TPZFMatrix< TVar >::Redim(), TPZMatrix< TVar >::Rows(), and UpdatePreviousState().
Referenced by SetAllCreateFunctionsContinuousWithMem(), and UpdatePreviousState().
|
overridevirtual |
Save the element data to a stream.
Save the element data to a stream
Reimplemented from TPZSavable.
Reimplemented in TPZSubCompMesh, TPZFlowCompMesh, and TPZCompMeshReferred.
Definition at line 1975 of file pzcmesh.cpp.
References fBlock, fConnectVec, fCreate, fDefaultOrder, fDimModel, fElementSolution, fElementVec, fGMesh, fMaterialVec, fName, fNmeshes, fReference, fSolution, fSolutionBlock, TPZStream::Write(), TPZCreateApproximationSpace::Write(), TPZAdmChunkVector< T, EXP >::Write(), TPZBlock< TVar >::Write(), TPZFMatrix< TVar >::Write(), TPZPersistenceManager::WritePointer(), and TPZStream::WritePointers().
Referenced by main(), SetAllCreateFunctionsContinuousWithMem(), TPZCompMeshReferred::Write(), TPZFlowCompMesh::Write(), and TPZSubCompMesh::Write().
|
protected |
Block structure to right construction of the stiffness matrix and load vector.
Definition at line 80 of file pzcmesh.h.
Referenced by AllocateNewConnect(), TPZSubCompMesh::AllocateNewConnect(), AutoBuildContDisc(), BandWidth(), Block(), CleanUp(), CleanUpUnconnectedNodes(), ExpandSolution(), TPZFlowCompMesh::ExpandSolution2(), TPZSubCompMesh::LoadSolution(), NEquations(), TPZSubCompMesh::NumberRigidBodyModes(), operator=(), Permute(), Read(), TPZSubCompMesh::SetNumberRigidBodyModes(), Skyline(), TPZCompMesh(), and Write().
|
protected |
List of pointers to nodes.
Definition at line 62 of file pzcmesh.h.
Referenced by AllocateNewConnect(), TPZSubCompMesh::AllocateNewConnect(), BandWidth(), TPZSubCompMesh::CalcStiff(), TPZMultiphysicsCompMesh::CleanElementsConnects(), CleanUp(), CleanUpUnconnectedNodes(), ComputeElGraph(), ComputeFillIn(), TPZSubCompMesh::ComputeNodElCon(), ComputeNodElCon(), TPZSubCompMesh::ComputePermutationInternalFirst(), ConnectVec(), TPZSubCompMesh::LoadSolution(), TPZSubCompMesh::MakeExternal(), TPZSubCompMesh::MakeInternal(), NEquations(), NIndependentConnects(), TPZSubCompMesh::NumberRigidBodyModes(), TPZSubCompMesh::NumInternalEquations(), operator=(), Permute(), TPZSubCompMesh::PermuteExternalConnects(), TPZSubCompMesh::PotentialInternal(), TPZSubCompMesh::Print(), Print(), Read(), TPZSubCompMesh::SetNumberRigidBodyModes(), Skyline(), Write(), and TPZSubCompMesh::~TPZSubCompMesh().
|
protected |
The object which defines the type of space being created.
Definition at line 92 of file pzcmesh.h.
Referenced by ApproxSpace(), AutoBuild(), AutoBuildContDisc(), Coarsen(), Discontinuous2Continuous(), operator=(), Read(), and Write().
|
protected |
Default order for all elements of this mesh.
Definition at line 89 of file pzcmesh.h.
Referenced by GetDefaultOrder(), operator=(), Read(), TPZCompMesh(), and Write().
|
protected |
Definition at line 86 of file pzcmesh.h.
Referenced by Dimension(), operator=(), Read(), SetReference(), TPZCompMesh(), and Write().
|
protected |
Solution vectors organized by element.
Definition at line 83 of file pzcmesh.h.
Referenced by ElementSolution(), operator=(), Read(), SetElementSolution(), and Write().
|
protected |
List of pointers to elements.
Definition at line 60 of file pzcmesh.h.
Referenced by AdjustBoundaryElements(), AssembleError(), BandWidth(), BuildTransferMatrix(), TPZMultiphysicsCompMesh::CleanElementsConnects(), CleanUp(), Coarsen(), CompareMesh(), ComputeElGraph(), ComputeFillIn(), ComputeNodElCon(), TPZFlowCompMesh::ComputeTimeStep(), TPZSubCompMesh::CreateGraphicalElement(), Discontinuous2Continuous(), Divide(), ElementVec(), EvaluateError(), GetRefPatches(), LoadReferences(), TPZCompMeshReferred::LoadReferred(), LoadSolution(), operator=(), Print(), Read(), Skyline(), TPZCompMesh(), TPZSubCompMesh::TransferElementFrom(), Write(), and TPZSubCompMesh::~TPZSubCompMesh().
|
protected |
Autopointer to the geometric mesh used in case the user has passed an autopointer.
Definition at line 54 of file pzcmesh.h.
Referenced by Read(), SetReference(), and Write().
|
protected |
Map of pointers to materials.
Definition at line 65 of file pzcmesh.h.
Referenced by BuildTransferMatrix(), BuildTransferMatrixDesc(), CleanUp(), TPZFlowCompMesh::CollectFluidMaterials(), ComputeElGraph(), CopyMaterials(), FindMaterial(), InsertMaterialObject(), MaterialVec(), Print(), ProjectSolution(), Read(), TPZSubCompMesh::Write(), and Write().
|
protected |
Grid name for model identification.
Definition at line 57 of file pzcmesh.h.
Referenced by Name(), operator=(), Print(), Read(), SetName(), TPZCompMesh(), and Write().
|
protected |
Definition at line 98 of file pzcmesh.h.
Referenced by GetNMeshes(), Read(), TPZCompMesh(), and Write().
|
protected |
Geometric grid to which this grid refers.
Definition at line 51 of file pzcmesh.h.
Referenced by operator=(), Read(), Reference(), SetReference(), TPZCompMesh(), and Write().
|
protected |
Solution at previous state.
Definition at line 76 of file pzcmesh.h.
Referenced by SolutionN(), and UpdatePreviousState().
|
protected |
Solution vector.
Definition at line 72 of file pzcmesh.h.
Referenced by CleanUp(), ExpandSolution(), TPZFlowCompMesh::ExpandSolution2(), TPZSubCompMesh::LoadSolution(), LoadSolution(), operator=(), Permute(), Read(), Solution(), TPZCompMesh(), TPZSubCompMesh::TransferMultiphysicsElementSolution(), UpdatePreviousState(), and Write().
|
protected |
Block structure of the solution vector ????
Definition at line 68 of file pzcmesh.h.
Referenced by CleanUp(), ExpandSolution(), TPZFlowCompMesh::ExpandSolution2(), operator=(), Permute(), Read(), TPZCompMesh(), and Write().