NeoPZ
pzcmesh.h
Go to the documentation of this file.
1 
6 #ifndef PZCMESHHPP
7 #define PZCMESHHPP
8 
9 
10 #include <stddef.h> // for NULL
11 #include <iostream> // for operator<<, string, cout, ostream
12 #include <map> // for map
13 #include <set> // for set
14 #include "pzadmchunk.h" // for TPZAdmChunkVector
15 #include "pzblock.h" // for TPZBlock
16 #include "TPZChunkVector.h" // for TPZChunkVector
17 #include "pzconnect.h" // for TPZConnect
18 #include "pzcreateapproxspace.h" // for TPZCreateApproximationSpace
19 #include "pzgmesh.h" // for TPZGeoMesh
20 #include "pzmatrix.h" // for TPZFMatrix, TPZMatrix
21 #include "pzreal.h" // for STATE, REAL
22 #include "TPZSavable.h" // for TPZSavable
23 #include "pzstack.h" // for TPZStack
24 #include "pzvec.h" // for TPZVec
25 #include "tpzautopointer.h" // for TPZAutoPointer
26 #include "pzcheckgeom.h" // for TPZCheckGeom
27 #include <functional>
28 #include "pzcompel.h"
29 
30 class TPZGeoEl;
31 class TPZMaterial;
32 class TPZStream;
33 template <class TVar> class TPZTransfer;
34 
35 
47 class TPZCompMesh : public virtual TPZSavable {
48 
49 protected:
52 
55 
57  std::string fName;
58 
63 
65  std::map<int, TPZMaterial * > fMaterialVec;
67  //TPZBlock<REAL> fSolutionBlock;
69 
71  //TPZFMatrix<REAL> fSolution;
73 
77 
79  //TPZBlock<REAL> fBlock;
81 
84 
85  /* @brief set the dimension of the simulation or the model */
86  int fDimModel;
87 
90 
93 
95  void ModifyPermute(TPZVec<int64_t> &permute, int64_t lagrangeq, int64_t maxeq);
96 
97  //quantity of meshes associated with a mesh multiphysics
98  int64_t fNmeshes;
99 
100 public:
101 
106  TPZCompMesh(TPZGeoMesh* gr=0);
110  TPZCompMesh(const TPZCompMesh &copy);
112  TPZCompMesh &operator=(const TPZCompMesh &copy);
113 
115  virtual ~TPZCompMesh();
116 
123  virtual REAL CompareMesh(int var, char *matname);
128  void GetRefPatches(std::set<TPZGeoEl *> &grpatch);
130  void GetNodeToElGraph(TPZVec<int64_t> &nodtoelgraph, TPZVec<int64_t> &nodtoelgraphinde,TPZStack<int64_t> &elgraph, TPZVec<int64_t> &elgraphindexx);
131 
133  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);
134 
136  void SetName(const std::string &nm);
137 
139  void SetDimModel(int dim) { fDimModel = dim; }
140 
142  void SetReference(TPZGeoMesh * gmesh);
143 
146 
148  int Dimension() const {return fDimModel;}
149 
151  std::string &Name() {return fName;}
152 
154  void CleanUp();
155 
163  int64_t NConnects() const {return fConnectVec.NElements();}
164 
166  int64_t NIndependentConnects();
167 
169  int64_t NElements() const {return fElementVec.NElements();}
170 
172  size_t NMaterials() const {return fMaterialVec.size();}
173 
175 
184 
185  TPZCompEl * Element(int64_t iel)
186  {
187  return fElementVec[iel];
188  }
189 
190  const TPZCompEl * Element(int64_t iel) const {
191  return fElementVec[iel];
192  }
193 
196 
199 
201 
203  std::map<int ,TPZMaterial * > &MaterialVec() { return fMaterialVec; }
204 
206  std::map<int ,TPZMaterial * > MaterialVec() const { return fMaterialVec; }
207 
209  TPZGeoMesh *Reference() const { return fReference; }
210 
212  //const TPZBlock<REAL> &Block() const { return fBlock;}
213  const TPZBlock<STATE> &Block() const { return fBlock;}
214 
217 
220 
223 
226 
241  virtual int64_t AllocateNewConnect(int nshape, int nstate, int order);
242 
247  virtual int64_t AllocateNewConnect(const TPZConnect &connect);
248 
254 
256  void InitializeBlock();
257 
259  virtual void ExpandSolution();
260 
270  void SetElementSolution(int64_t i, TPZVec<STATE> &sol);
271 
284  virtual void Print(std::ostream & out = std::cout) const;
285 
290  void ConnectSolution(std::ostream & out);
291 
298  TPZMaterial * FindMaterial(int id);
299 
301  void LoadReferences();
302 
304  virtual void ComputeNodElCon();
305 
307  virtual void ComputeNodElCon(TPZVec<int> &nelconnected) const;
308 
310  virtual void CleanUpUnconnectedNodes();
311 
317  void ComputeElGraph(TPZStack<int64_t> &elgraph, TPZVec<int64_t> &elgraphindex);
318 
324  void ComputeElGraph(TPZStack<int64_t> &elgraph, TPZVec<int64_t> &elgraphindex, std::set<int> & mat_ids);
325 
327  virtual void BuildCornerConnectList(std::set<int64_t> &connectindexes) const;
328 
329 
337  virtual TPZCompMesh *FatherMesh() const { return NULL; }
338 
343  virtual void MakeInternal(int64_t local) {;}
344 
346  virtual void MakeAllInternal() {;}
347 
349  virtual TPZCompMesh* RootMesh (int64_t local) { return this; }
350 
356  virtual int64_t TransferElementFrom(TPZCompMesh *mesh, int64_t elindex) { return elindex; }
357 
363  virtual int64_t TransferElementTo(TPZCompMesh *mesh, int64_t elindex) { return elindex; }
364 
370  virtual int64_t TransferElement(TPZCompMesh *mesh, int64_t elindex) { return elindex; }
371 
378  virtual int64_t PutinSuperMesh (int64_t local, TPZCompMesh *super);
379 
386  virtual int64_t GetFromSuperMesh (int64_t superind, TPZCompMesh *super);
387 
392  virtual TPZCompMesh * CommonMesh (TPZCompMesh *mesh);
393 
394 
395 
399  {
400  return fDefaultOrder;
401  }
402 
403  void SetDefaultOrder( int order )
404  {
405  fDefaultOrder = order;
406  }
407 
415  int64_t NEquations();
416 
418  int BandWidth();
419 
424  virtual void Skyline(TPZVec<int64_t> &skyline);
425 
431  void AssembleError(TPZFMatrix<REAL> &estimator, int errorid);
432 
438  void BuildTransferMatrix(TPZCompMesh &coarsemesh, TPZTransfer<STATE> &transfer);
439 
441  void BuildTransferMatrixDesc(TPZCompMesh &transfermesh,TPZTransfer<STATE> &transfer);
442  void ProjectSolution(TPZFMatrix<STATE> &projectsol);
443 
444  //set nummber of meshs
445  void SetNMeshes(int64_t nmeshes){
446  fNmeshes = nmeshes;
447  }
448 
449  int64_t GetNMeshes(){
450  return fNmeshes;
451  }
452 
453 private:
454 
457  virtual void AutoBuild(const std::set<int> *MaterialIDs);
458 
459 public:
460 
462  TPZCompEl *CreateCompEl(TPZGeoEl *gel, int64_t &index)
463  {
464  return fCreate.CreateCompEl(gel, *this, index);
465  }
466 
469  virtual void AutoBuild(const std::set<int> &MaterialIDs){
470  fCreate.BuildMesh(*this,MaterialIDs);
471  }
472 
474  virtual void AutoBuild(){
475 #ifdef PZDEBUG
476  {
477  TPZGeoMesh *gmesh = Reference();
478  TPZCheckGeom check(gmesh);
479  check.CheckUniqueId();
480  }
481 #endif
482  fCreate.BuildMesh(*this);
483  }
484 
486  void AutoBuild(const TPZVec<int64_t> &gelindexes)
487  {
488  fCreate.BuildMesh(*this,gelindexes);
489  }
490 
496  virtual void AutoBuildContDisc(const TPZVec<TPZGeoEl*> &continuous, const TPZVec<TPZGeoEl*> &discontinuous);
497 
499  {
500  return fCreate;
501  }
502 
504  {
506  }
507 
509  {
511  }
512 
514  {
516  }
517 
519  {
521  }
522 
524  {
526  }
527 
528 
529 
530 #ifndef STATE_COMPLEX
532  {
534  }
535 #endif
536 
538  {
539  fCreate.SetAllCreateFunctions(cel, this);
540  }
541 
543  {
545  }
546 
548  {
550  }
551 
553  {
555  }
556 
559  int Consolidate();
560 
564  int Check();
565 
570  void LoadSolution(const TPZFMatrix<STATE> &sol);
571 
574  void UpdatePreviousState(STATE mult = 1.0);
575 
580 
587  void Divide(int64_t index,TPZVec<int64_t> &subindex,int interpolate = 0);
588 
595  void Coarsen(TPZVec<int64_t> &elements, int64_t &index, bool CreateDiscontinuous = false);
596 
599 
604  void AdjustBoundaryElements();
605 
611  void Permute(TPZVec<int64_t> &permute);
612 
614  void SaddlePermute();
615  void SaddlePermute2();
616 
618  /*
619  * @param varname name of the variable that will be integrated
620  * @param matids ids of the materials that will contribute to the integral
621  */
622  TPZVec<STATE> Integrate(const std::string &varname, const std::set<int> &matids);
623 
624 
633  void EvaluateError(std::function<void (const TPZVec<REAL> &loc,TPZVec<STATE> &val,TPZFMatrix<STATE> &deriv)> fp, bool store_error,
634  TPZVec<REAL> &errorSum);
635 
647  void ConvertDiscontinuous2Continuous(REAL eps, int opt, int dim, TPZVec<STATE> &celJumps);
648 
654  void Discontinuous2Continuous(int64_t disc_index, int64_t &new_index);
655 
660  TPZCompMesh * Clone() const;
661 
663  void CopyMaterials(TPZCompMesh &mesh) const ;
664 
665  REAL DeltaX();
666 
667  REAL MaximumRadiusOfMesh();
668 
669  REAL LesserEdgeOfMesh();
670 
672  TPZCompMesh *ComputeMesh(TPZVec<int64_t> &accumlist,int64_t numaggl);
673 
679  void ComputeFillIn(int64_t resolution, TPZFMatrix<REAL> &fillin);
680 
682  public:
683 int ClassId() const override;
684 
686  void Write(TPZStream &buf, int withclassid) const override;
687 
689  void Read(TPZStream &buf, void *context) override;
690 
691 };
692 
693 inline int64_t TPZCompMesh::AllocateNewConnect(int nshape, int nstate, int order) {
694  int64_t connectindex = fConnectVec.AllocateNewElement();
695  TPZConnect &c = fConnectVec[connectindex];
696  c.SetNShape(nshape);
697  c.SetNState(nstate);
698  c.SetOrder(order,connectindex);
700  int64_t blocknum = fBlock.NBlocks();
701  fBlock.SetNBlocks(blocknum+1);
702  fBlock.Set(blocknum,nshape*nstate);
703  c.SetSequenceNumber(blocknum);
704  return connectindex;
705 }
706 
711 inline int64_t TPZCompMesh::AllocateNewConnect(const TPZConnect &connect)
712 {
713  int64_t connectindex = fConnectVec.AllocateNewElement();
714  TPZConnect &c = fConnectVec[connectindex];
715  c = connect;
716  c.RemoveDepend();
717  int nshape = c.NShape();
718  int nstate = c.NState();
719  int64_t blocknum = fBlock.NBlocks();
720  fBlock.SetNBlocks(blocknum+1);
721  fBlock.Set(blocknum,nshape*nstate);
722  c.SetSequenceNumber(blocknum);
723  return connectindex;
724 }
725 
726 
728  this->fReference = gmesh;
729  this->fDimModel = gmesh->Dimension();
731 }
732 
734  this->fReference = gmesh.operator->();
735  this->fDimModel = gmesh->Dimension();
736  fGMesh = gmesh;
737 }
738 
739 #endif
int64_t NElements() const
Number of computational elements allocated.
Definition: pzcmesh.h:169
virtual void MakeInternal(int64_t local)
Makes a specified connection a internal mesh connection.
Definition: pzcmesh.h:343
void SetAllCreateFunctionsHDivPressure(int meshdim)
Create an approximation space with HDiv elements and full basis for quadrilateral element...
TPZFMatrix< STATE > & ElementSolution()
Access method for the element solution vectors.
Definition: pzcmesh.h:225
TPZGeoMesh * Reference() const
Returns a pointer to the geometrical mesh associated.
Definition: pzcmesh.h:209
void SetAllCreateFunctionsMultiphysicElemWithMem()
Create an approximation space based on multiphysics elements with memory.
void SetAllCreateFunctionsDiscontinuous()
Create discontinuous approximation spaces.
This class performs a series of consistency tests on geometric transformations between elements...
Definition: pzcheckgeom.h:16
void SetAllCreateFunctions(TPZCompEl &cel, TPZCompMesh *mesh)
Create approximation spaces corresponding to the space defined by cel.
int fDimModel
Definition: pzcmesh.h:86
void SetElementSolution(int64_t i, TPZVec< STATE > &sol)
Set a ith element solution, expanding the element-solution matrix if necessary.
Definition: pzcmesh.cpp:1533
int AllocateNewElement()
Makes more room for new elements.
Definition: pzadmchunk.h:184
virtual TPZCompMesh * FatherMesh() const
Get the father meshes stack.
Definition: pzcmesh.h:337
void AdjustBoundaryElements()
Will refine the elements associated with a boundary condition till there are no elements constrained ...
Definition: pzcmesh.cpp:1423
Contains declaration of the TPZSavable class which defines the interface to save and restore objects ...
void SetAllCreateFunctionsContinuousWithMem()
Create an approximation space with continous elements with memory. Only dimension 3 elements quem hav...
void SetAllCreateFunctionsContinuous()
Definition: pzcmesh.h:508
int Set(const int index, const int dim, const int pos=-1)
Modifies existing block dimensions or creates a new block with given index.
Definition: pzblock.cpp:104
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
void SetName(const std::string &nm)
Set the mesh name.
Definition: pzcmesh.cpp:232
TPZFMatrix< STATE > fSolution
Solution vector.
Definition: pzcmesh.h:72
TPZAdmChunkVector< TPZCompEl * > fElementVec
List of pointers to elements.
Definition: pzcmesh.h:60
void SetOrder(int order, int64_t index)
Set the order of the shapefunction associated with the connect.
Definition: pzconnect.h:168
int64_t NEquations()
This computes the number of equations associated with non-restrained nodes.
Definition: pzcmesh.cpp:721
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 ...
Definition: pzcmesh.cpp:2052
int64_t NIndependentConnects()
Number of independent connect objects.
Definition: pzcmesh.cpp:1080
void CheckUniqueId()
Verify is the ids of the elements and nodes are unique.
virtual int64_t AllocateNewConnect(int nshape, int nstate, int order)
Returns an index to a new connect.
Definition: pzcmesh.h:693
void SaddlePermute()
Put the sequence number of the pressure connects after the seq number of the flux connects...
Definition: pzcmesh.cpp:2328
int64_t fNmeshes
Definition: pzcmesh.h:98
Contains declaration of TPZCompEl class which defines the interface of a computational element...
Templated vector implementation.
TPZCompMesh * ComputeMesh(TPZVec< int64_t > &accumlist, int64_t numaggl)
Creates a mesh inserting into accumlist the element list to more refined mesh.
int ClassId() const override
Returns the unique identifier for reading/writing objects to streams.
Definition: pzcmesh.cpp:1969
Contains declaration of TPZCheckGeom class which performs a series of consistency tests on geometric ...
const TPZAdmChunkVector< TPZConnect > & ConnectVec() const
Definition: pzcmesh.h:200
TPZVec< STATE > Integrate(const std::string &varname, const std::set< int > &matids)
Integrate the variable name over the mesh.
Definition: pzcmesh.cpp:2162
std::string & Name()
Returns the mesh name.
Definition: pzcmesh.h:151
TPZFMatrix< STATE > & SolutionN()
Access the previous solution vector.
Definition: pzcmesh.h:222
void SetLagrangeMultiplier(unsigned char mult)
Set the connect as a pressure connect or not.
Definition: pzconnect.h:236
REAL LesserEdgeOfMesh()
Definition: pzcmesh.cpp:1849
std::map< int,TPZMaterial *> MaterialVec() const
Returns a reference to the material pointers vector.
Definition: pzcmesh.h:206
void SetAllCreateFunctionsContinuousReferred()
Create a continuous approximation space with referred elements.
Declarates the TPZBlock<REAL>class which implements block matrices.
virtual int64_t TransferElement(TPZCompMesh *mesh, int64_t elindex)
Transfer one element form a submesh to another mesh.
Definition: pzcmesh.h:370
size_t NMaterials() const
Number of materials.
Definition: pzcmesh.h:172
int64_t NElements() const
Access method to query the number of elements of the vector.
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
virtual int64_t TransferElementFrom(TPZCompMesh *mesh, int64_t elindex)
Transfer one element from a specified mesh to the current submesh.
Definition: pzcmesh.h:356
Implements a chunk vector with free store administration. Utility.
Definition: TPZStream.h:39
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 globa...
Definition: pzcmesh.cpp:1869
virtual void BuildCornerConnectList(std::set< int64_t > &connectindexes) const
adds the connect indexes associated with base shape functions to the set
Definition: pzcmesh.cpp:2782
void CleanUp()
Delete all the dynamically allocated data structures.
Definition: pzcmesh.cpp:161
void SetAllCreateFunctionsDiscontinuousReferred()
Create a discontinuous approximation space with referred elements.
void LoadReferences()
Map this grid in the geometric grid.
Definition: pzcmesh.cpp:482
void RemakeAllInterfaceElements()
Deletes all interfaces and rebuild them all.
Definition: pzcmesh.cpp:1281
void ComputeElGraph(TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindex)
Computes the connectivity graph of the elements, as appropriate for the TPZRenumbering class...
Definition: pzcmesh.cpp:1091
void SetAllCreateFunctionsMultiphysicElem()
Create an approximation space based on multiphysics elements.
REAL DeltaX()
Definition: pzcmesh.cpp:1810
void SetDefaultOrder(int order)
Definition: pzcmesh.h:403
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
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.
Definition: pzcmesh.cpp:1599
void BuildTransferMatrix(TPZCompMesh &coarsemesh, TPZTransfer< STATE > &transfer)
Builds the transfer matrix from the current grid to the coarse grid.
Definition: pzcmesh.cpp:883
TPZAutoPointer< TPZGeoMesh > fGMesh
Autopointer to the geometric mesh used in case the user has passed an autopointer.
Definition: pzcmesh.h:54
void SetAllCreateFunctionsHDiv(int meshdim)
Create an approximation space with HDiv elements.
virtual void AutoBuildContDisc(const TPZVec< TPZGeoEl *> &continuous, const TPZVec< TPZGeoEl *> &discontinuous)
Creates the computational elements, and the degree of freedom nodes.
Definition: pzcmesh.cpp:326
TPZCompEl * CreateCompEl(TPZGeoEl *gel, int64_t &index)
Create a computational element based on the geometric element.
Definition: pzcmesh.h:462
TPZCompEl * Element(int64_t iel)
Definition: pzcmesh.h:185
void LoadSolution(const TPZFMatrix< STATE > &sol)
Given the solution of the global system of equations, computes and stores the solution for the restri...
Definition: pzcmesh.cpp:441
int fDefaultOrder
Default order for all elements of this mesh.
Definition: pzcmesh.h:89
TPZCreateApproximationSpace & ApproxSpace()
Definition: pzcmesh.h:498
Contains declaration of TPZMesh class which defines a geometrical mesh and contains a corresponding l...
void GetRefPatches(std::set< TPZGeoEl *> &grpatch)
Gives all patches of the mesh.
Definition: pzcmesh.cpp:1558
void SetAllCreateFunctionsDiscontinuousReferred()
Definition: pzcmesh.h:513
virtual void Skyline(TPZVec< int64_t > &skyline)
This method computes the skyline of the system of equations.
Definition: pzcmesh.cpp:787
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
const TPZAdmChunkVector< TPZCompEl * > & ElementVec() const
Returns a reference to the element pointers vector.
Definition: pzcmesh.h:195
TPZCompMesh * Clone() const
Clone this mesh.
Definition: pzcmesh.cpp:1720
void SetAllCreateFunctionsDiscontinuous()
Definition: pzcmesh.h:503
void Read(TPZStream &buf, void *context) override
Read the element data from a stream.
Definition: pzcmesh.cpp:2006
void SetAllCreateFunctions(TPZCompEl &cel)
Definition: pzcmesh.h:537
virtual void ExpandSolution()
Adapt the solution vector to new block dimensions.
Definition: pzcmesh.cpp:396
void SetSequenceNumber(int64_t i)
Set the sequence number for the global system of equations of the connect object. ...
Definition: pzconnect.h:165
void SetAllCreateFunctionsContinuousWithMem()
Definition: pzcmesh.h:552
void SetAllCreateFunctionsHDiv()
Definition: pzcmesh.h:523
Contains declaration of TPZConnect class which represents a set of shape functions associated with a ...
void SetNState(int nstate)
Set the number of state variables.
Definition: pzconnect.h:196
int Dimension() const
Returns the dimension of the simulation.
Definition: pzcmesh.h:148
void Coarsen(TPZVec< int64_t > &elements, int64_t &index, bool CreateDiscontinuous=false)
Create a computational element father for the comp. elem. into elements.
Definition: pzcmesh.cpp:1164
void SetReference(TPZGeoMesh *gmesh)
Sets the geometric reference mesh.
Definition: pzcmesh.h:727
virtual void ComputeNodElCon()
Compute the number of elements connected to each connect object.
Definition: pzcmesh.cpp:668
void Divide(int64_t index, TPZVec< int64_t > &subindex, int interpolate=0)
Divide the element corresponding to index.
Definition: pzcmesh.cpp:1146
void SaddlePermute2()
Definition: pzcmesh.cpp:2588
TPZMaterial * FindMaterial(int id)
Find the material with identity id.
Definition: pzcmesh.cpp:297
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...
Definition: pzcmesh.cpp:1512
unsigned int NShape() const
Definition: pzconnect.h:151
void SetNShape(int nshape)
Set the number of shape functions associated with the connect.
Definition: pzconnect.h:205
void GetNodeToElGraph(TPZVec< int64_t > &nodtoelgraph, TPZVec< int64_t > &nodtoelgraphinde, TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindexx)
Gives the conects graphs.
Definition: pzcmesh.cpp:1576
TPZBlock< STATE > & Block()
Access the block structure of the solution vector.
Definition: pzcmesh.h:216
TPZFMatrix< STATE > fElementSolution
Solution vectors organized by element.
Definition: pzcmesh.h:83
void SetAllCreateFunctionsContinuous()
Create continuous approximation spaces.
virtual int64_t TransferElementTo(TPZCompMesh *mesh, int64_t elindex)
Transfer one element from a submesh to another mesh.
Definition: pzcmesh.h:363
const TPZBlock< STATE > & Block() const
Access the block structure of the solution vector.
Definition: pzcmesh.h:213
Contains declaration of the TPZAutoPointer class which has Increment and Decrement actions are mutexe...
void SetNMeshes(int64_t nmeshes)
Definition: pzcmesh.h:445
int64_t NConnects() const
Number of connects allocated including free nodes.
Definition: pzcmesh.h:163
TPZCompMesh(TPZGeoMesh *gr=0)
Constructor from geometrical mesh.
Definition: pzcmesh.cpp:63
TPZGeoMesh * fReference
Geometric grid to which this grid refers.
Definition: pzcmesh.h:51
TPZAdmChunkVector< TPZConnect > & ConnectVec()
Return a reference to the connect pointers vector.
Definition: pzcmesh.h:198
TPZFMatrix< STATE > fSolN
Solution at previous state.
Definition: pzcmesh.h:76
int SetNBlocks(const int num_of_blocks)
Sets number of blocks on diagonal matrix.
Definition: pzblock.cpp:91
int Consolidate()
Will build the list of element boundary conditions build the list of connect boundary conditions...
TPZBlock< STATE > fBlock
Block structure to right construction of the stiffness matrix and load vector.
Definition: pzcmesh.h:80
int BandWidth()
This method computes the bandwidth of the system of equations.
Definition: pzcmesh.cpp:747
void SetDimModel(int dim)
Set de dimension of the domain of the problem.
Definition: pzcmesh.h:139
unsigned char NState() const
Number of state variables associated with the connect.
Definition: pzconnect.h:146
void UpdatePreviousState(STATE mult=1.0)
Definition: pzcmesh.cpp:2826
A simple stack.
std::map< int,TPZMaterial *> & MaterialVec()
Returns a reference to the material pointers vector.
Definition: pzcmesh.h:203
Contains TPZMatrix<TVar>class, root matrix class.
int NBlocks() const
Returns number of blocks on diagonal.
Definition: pzblock.h:165
TPZBlock< STATE > fSolutionBlock
Block structure of the solution vector ????
Definition: pzcmesh.h:68
virtual void MakeAllInternal()
Make all mesh connections internal mesh connections. Connects to an internal connection.
Definition: pzcmesh.h:346
int InsertMaterialObject(TPZMaterial *mat)
Insert a material object in the datastructure.
Definition: pzcmesh.cpp:287
virtual void CleanUpUnconnectedNodes()
Delete the nodes which have no elements connected to them.
Definition: pzcmesh.cpp:498
TPZCompEl * CreateCompEl(TPZGeoEl *gel, TPZCompMesh &mesh, int64_t &index) const
Create a computational element using the function pointer for the topology.
void AssembleError(TPZFMatrix< REAL > &estimator, int errorid)
Assemble the vector with errors estimators.
Definition: pzcmesh.cpp:2113
virtual ~TPZCompMesh()
Simple Destructor.
Definition: pzcmesh.cpp:134
int Dimension()
Get Dimension.
Definition: pzgmesh.h:291
Free store vector implementation in chunks.
This class implements a geometric mesh for the pz environment. Geometry.
Definition: pzgmesh.h:48
virtual REAL CompareMesh(int var, char *matname)
This method will initiate the comparison between the current computational mesh and the mesh which is...
Definition: pzcmesh.cpp:1522
const TPZCompEl * Element(int64_t iel) const
Definition: pzcmesh.h:190
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
void SetAllCreateFunctionsHDivPressure()
Definition: pzcmesh.h:531
Contains the functions to create different computational elements (one- two- three-dimensional).
TPZAdmChunkVector< TPZCompEl * > & ElementVec()
Returns a reference to the element pointers vector.
Definition: pzcmesh.h:183
void BuildTransferMatrixDesc(TPZCompMesh &transfermesh, TPZTransfer< STATE > &transfer)
To discontinuous elements.
Definition: pzcmesh.cpp:941
void BuildMesh(TPZCompMesh &cmesh, const std::set< int > &MaterialIDs) const
Creates the computational elements, and the degree of freedom nodes.
virtual TPZCompMesh * CommonMesh(TPZCompMesh *mesh)
Gives the commom father mesh of the specified mesh and the current submesh.
Definition: pzcmesh.cpp:2794
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.
Definition: pzcmesh.cpp:1395
TPZCompMesh & operator=(const TPZCompMesh &copy)
copy the content of the mesh
Definition: pzcmesh.cpp:1683
TPZCreateApproximationSpace fCreate
The object which defines the type of space being created.
Definition: pzcmesh.h:92
virtual void AutoBuild(const std::set< int > &MaterialIDs)
Creates the computational elements, and the degree of freedom nodes.
Definition: pzcmesh.h:469
std::string fName
Grid name for model identification.
Definition: pzcmesh.h:57
void RemoveDepend(int64_t myindex, int64_t dependindex)
Remove dependency between connects if exist.
Definition: pzconnect.cpp:178
void CopyMaterials(TPZCompMesh &mesh) const
Copies the materials of this mesh to the given mesh.
Definition: pzcmesh.cpp:1797
void AutoBuild(const TPZVec< int64_t > &gelindexes)
build the computational elements for the geometric element indexes
Definition: pzcmesh.h:486
std::map< int, TPZMaterial *> fMaterialVec
Map of pointers to materials.
Definition: pzcmesh.h:65
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
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...
Definition: pzcmesh.cpp:2733
void SetAllCreateFunctionsContinuousReferred()
Definition: pzcmesh.h:518
void SetAllCreateFunctionsMultiphysicElemWithMem()
Definition: pzcmesh.h:547
virtual TPZCompMesh * RootMesh(int64_t local)
Returns the rootmesh who have the specified connection.
Definition: pzcmesh.h:349
void InitializeBlock()
Resequence the block object, remove unconnected connect objects and reset the dimension of the soluti...
Definition: pzcmesh.cpp:391
Implements rectangular matrix which extends a solution vector of the coarse mesh to a solution vector...
Definition: pzcmesh.h:33
void ProjectSolution(TPZFMatrix< STATE > &projectsol)
Definition: pzcmesh.cpp:1923
REAL MaximumRadiusOfMesh()
Definition: pzcmesh.cpp:1833
This class defines the interface to save and restore objects from TPZStream objects. Persistency.
Definition: TPZSavable.h:67
void ConnectSolution(std::ostream &out)
Print the solution by connect index.
Definition: pzcmesh.cpp:1373
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
void SetAllCreateFunctionsMultiphysicElem()
Definition: pzcmesh.h:542
int GetDefaultOrder()
Definition: pzcmesh.h:398
void TransferMultiphysicsSolution()
Transfer multiphysics mesh solution.
Definition: pzcmesh.cpp:470
void Discontinuous2Continuous(int64_t disc_index, int64_t &new_index)
This method convert a discontinuous element with index disc_index in continuous element.
Definition: pzcmesh.cpp:1234
TPZFMatrix< STATE > & Solution()
Access the solution vector.
Definition: pzcmesh.h:219
virtual void AutoBuild()
Creates the computational elements, and the degree of freedom nodes.
Definition: pzcmesh.h:474
void Write(TPZStream &buf, int withclassid) const override
Save the element data to a stream.
Definition: pzcmesh.cpp:1975
TPZAdmChunkVector< TPZConnect > fConnectVec
List of pointers to nodes.
Definition: pzcmesh.h:62
int64_t GetNMeshes()
Definition: pzcmesh.h:449
virtual void Print(std::ostream &out=std::cout) const
Prints mesh data.
Definition: pzcmesh.cpp:236
virtual int64_t GetFromSuperMesh(int64_t superind, TPZCompMesh *super)
Get an external connection from the supermesh - Supermesh is one mesh who contains the analised subme...
Definition: pzcmesh.cpp:1517
void Permute(TPZVec< int64_t > &permute)
Permute the sequence number of the connect objects It is a permute gather operation.
Definition: pzcmesh.cpp:1321