NeoPZ
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | List of all members
TPZDohrSubstruct< TVar > Class Template Reference

Implements sub structure matrices using Dohrman algorithm. Sub Structure. More...

#include <tpzdohrsubstruct.h>

Inheritance diagram for TPZDohrSubstruct< TVar >:
[legend]
Collaboration diagram for TPZDohrSubstruct< TVar >:
[legend]

Public Types

enum  EWeightType { CorrectWeight, Uniform }
 

Public Member Functions

 TPZDohrSubstruct ()
 
 ~TPZDohrSubstruct ()
 
int ClassId () const override
 Define the class id associated with the class. More...
 
void ReallocMatRed ()
 
void ContributeResidual (TPZFMatrix< TVar > &u, TPZFMatrix< TVar > &r)
 ContributeResidual. More...
 
void ContributeTestV1 (TPZFMatrix< TVar > &testV1, int NumCoarse)
 
void LoadWeightedResidual (const TPZFMatrix< TVar > &r_global)
 It computes W(i)*R(i)*r and stores it in fLocalWeightedResidual;. More...
 
void AdjustResidual (TPZFMatrix< TVar > &r_global)
 Adjust the residual to reflect a static condensation. More...
 
void Contribute_rc (TPZFMatrix< TVar > &rc)
 It computes the local contribution to r(c). More...
 
void Contribute_rc_local (TPZFMatrix< TVar > &residual_local, TPZFMatrix< TVar > &rc_local) const
 It computes the local contribution to r(c). More...
 
void Contribute_Kc (TPZMatrix< TVar > &Kc, TPZVec< int > &coarseindex)
 It computes the local contribution to K(c) More...
 
void Contribute_v1 (TPZFMatrix< TVar > &v1, TPZFMatrix< TVar > &invKc_rc)
 It computes the local contribution to v1. More...
 
void Contribute_v1_local (TPZFMatrix< TVar > &v1_local, TPZFMatrix< TVar > &invKc_rc) const
 It computes the local contribution to v1. More...
 
void Contribute_v2 (TPZFMatrix< TVar > &v2)
 It computes the local contribution to v2. More...
 
void Contribute_v2_local (TPZFMatrix< TVar > &residual_local, TPZFMatrix< TVar > &v2_local)
 It computes the local contribution to v2. More...
 
void Contribute_v3 (TPZFMatrix< TVar > &v3, const TPZFMatrix< TVar > &r, TPZFMatrix< TVar > &v1Plusv2) const
 It computes the local contribution to v(3) More...
 
void Contribute_v3_local (TPZFMatrix< TVar > &v3, TPZFMatrix< TVar > &v1Plusv2) const
 It computes the local contribution to v(3) More...
 
void Print (std::ostream &out) const
 
void ContributeKU (const TVar alpha, const TPZFMatrix< TVar > &uglobal, TPZFMatrix< TVar > &z) const
 Contribute to the global matrix vector multiplication. It's needed to the MultAdd method of TPZDohrMatrix. More...
 
void ContributeKULocal (const TVar alpha, const TPZFMatrix< TVar > &u, TPZFMatrix< TVar > &z) const
 Compute the multiplication of the local stiffness matrix with the vector u. More...
 
void ContributeGlobalDiagonal (TPZFMatrix< TVar > &StiffnessDiag)
 Computes the contribution of each substructure node to global Stiffness diagonal (or something like that). More...
 
void ContributeDiagonalLocal (TPZFMatrix< TVar > &StiffnessDiag)
 Computes the contribution of each substructure node to global Stiffness diagonal (or something like that). More...
 
void ComputeWeights (TPZFMatrix< TVar > &StiffnessDiag)
 Computes the weight matrix. More...
 
void ComputeWeightsLocal (TPZFMatrix< TVar > &StiffnessDiagLocal)
 Computes the weight matrix. More...
 
void Initialize ()
 Initializes the substructure. More...
 
void GetCoarseResidual (TPZFMatrix< TVar > &rc, TPZVec< int > &indices)
 Assembles the contribution to the coarse residual. More...
 
void GetCoarseStiffness (TPZMatrix< TVar > &stiff, TPZVec< int > &indices)
 Assembles the coarse dof stiffness matrix. More...
 
void PrepareSystems ()
 It prepares the datas for solving systems for phi and zi. More...
 
void SolveSystemPhi ()
 Solves the system for Phi and for v2. More...
 
void SolveSystemZi ()
 Solves the system for zi. More...
 
void ComputeCoarseStiffness ()
 Computes K(ci) and stores it in fKCi. More...
 
void AddInternalSolution (TPZFMatrix< TVar > &sol)
 Add the internal solution to the final result. More...
 
int NumEquations () const
 
- 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 void Write (TPZStream &buf, int withclassid) const
 Writes this object to the TPZStream buffer. Include the classid if withclassid = true. More...
 
virtual void Read (TPZStream &buf, void *context)
 read objects from the stream More...
 
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
 

Public Attributes

TPZFMatrix< TVar > fC_star
 Variables needed to solve the systems for phi and zi. More...
 
TPZFMatrix< TVar > fKeC_star
 
TPZStepSolver< TVar > finv
 
TPZFMatrix< TVar > fNullPivots
 
int fNEquations
 Number of equations of the substructure. More...
 
TPZVec< int > fCoarseNodes
 
TPZFMatrix< TVar > fC
 Constraint definition associated with this substructure. More...
 
TPZFMatrix< TVar > fPhiC
 Vectors associated with each constraint. More...
 
TPZFMatrix< TVar > fPhiC_Weighted_Condensed
 Phi * W matrix and condensed to the equations which are part of the global system. More...
 
TPZVec< int > fCoarseIndex
 Global index associated with each coarse degree of freedom. More...
 
TPZManVector< TVar > fWeights
 Weights associated with each variable/equation. More...
 
TPZAutoPointer< TPZMatrix< TVar > > fStiffness
 Stiffness matrix associated with the substructure. More...
 
TPZFMatrix< TVar > fKCi
 Stiffness matrix associated with the constraints. More...
 
TPZVec< int > fGlobalIndex
 Global vector indices of the equations/degrees of freedom. More...
 
TPZVec< std::pair< int, int > > fGlobalEqs
 Indices of the equations in the global system corresponding to the local equations. More...
 
TPZStepSolver< TVar > fInvertedStiffness
 Inverted (LU or Cholesky) stiffness matrix. More...
 
TPZStepSolver< TVar > fInvertedInternalStiffness
 Inverted (LU or Cholesky or LDLt) stiffness matrix for the internal degrees of freedom. More...
 
TPZStepSolver< TVar > fKCInvert
 Inverted restraint matrix. More...
 
TPZVec< int > fInternalEqs
 Internal nodes. More...
 
TPZVec< int > fBoundaryEqs
 Equations corresponding to boundary of the substructure. More...
 
TPZFMatrix< TVar > fLocalLoad
 Local load vector. More...
 
TPZFMatrix< TVar > fLocalWeightedResidual
 Local weighted residual - $ W(i)*R(i)*r $. More...
 
TPZFMatrix< TVar > fzi
 Needed to compute v2. Is the solution of a system. More...
 
TPZFMatrix< TVar > fAdjustSolution
 Solution vector which needs to be added to the converged system solution. More...
 

Static Public Attributes

static EWeightType fWeightType
 

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 TPZSavableCreateInstance (const int &classId)
 

Detailed Description

template<class TVar>
class TPZDohrSubstruct< TVar >

Implements sub structure matrices using Dohrman algorithm. Sub Structure.

Author
Philippe Devloo

Definition at line 28 of file tpzdohrsubstruct.h.

Member Enumeration Documentation

◆ EWeightType

template<class TVar>
enum TPZDohrSubstruct::EWeightType
Enumerator
CorrectWeight 
Uniform 

Definition at line 35 of file tpzdohrsubstruct.h.

Constructor & Destructor Documentation

◆ TPZDohrSubstruct()

template<class TVar >
TPZDohrSubstruct< TVar >::TPZDohrSubstruct ( )

Definition at line 21 of file tpzdohrsubstruct.cpp.

◆ ~TPZDohrSubstruct()

template<class TVar >
TPZDohrSubstruct< TVar >::~TPZDohrSubstruct ( )

Definition at line 27 of file tpzdohrsubstruct.cpp.

Member Function Documentation

◆ AddInternalSolution()

template<class TVar >
void TPZDohrSubstruct< TVar >::AddInternalSolution ( TPZFMatrix< TVar > &  sol)

Add the internal solution to the final result.

Add the internal solution to the final result

Definition at line 905 of file tpzdohrsubstruct.cpp.

References TPZMatrix< TVar >::Cols().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ AdjustResidual()

template<class TVar >
void TPZDohrSubstruct< TVar >::AdjustResidual ( TPZFMatrix< TVar > &  r_global)

Adjust the residual to reflect a static condensation.

The residual corresponding to the internal nodes will be zeroed

Adjust the residual to reflect a static condensation The residual corresponding to the internal nodes will be zeroed

Definition at line 847 of file tpzdohrsubstruct.cpp.

References fabs, and LOGPZ_DEBUG.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ClassId()

template<class TVar>
int TPZDohrSubstruct< TVar >::ClassId ( ) const
inlineoverridevirtual

Define the class id associated with the class.

This id has to be unique for all classes A non unique id is flagged at the startup of the program

Implements TPZSavable.

Definition at line 37 of file tpzdohrsubstruct.h.

References Hash().

◆ ComputeCoarseStiffness()

template<class TVar >
void TPZDohrSubstruct< TVar >::ComputeCoarseStiffness ( )

Computes K(ci) and stores it in fKCi.

Definition at line 430 of file tpzdohrsubstruct.cpp.

References LOGPZ_DEBUG, and TPZFMatrix< TVar >::Resize().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ComputeWeights()

template<class TVar >
void TPZDohrSubstruct< TVar >::ComputeWeights ( TPZFMatrix< TVar > &  StiffnessDiag)

Computes the weight matrix.

Parameters
StiffnessDiagis the diagonal of the global matrix (global numbering)

Definition at line 533 of file tpzdohrsubstruct.cpp.

References LOGPZ_DEBUG.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ComputeWeightsLocal()

template<class TVar >
void TPZDohrSubstruct< TVar >::ComputeWeightsLocal ( TPZFMatrix< TVar > &  StiffnessDiagLocal)

Computes the weight matrix.

Parameters
StiffnessDiagLocalis the diagonal of the global matrix (local numbering)

Definition at line 555 of file tpzdohrsubstruct.cpp.

References LOGPZ_DEBUG.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_Kc()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_Kc ( TPZMatrix< TVar > &  Kc,
TPZVec< int > &  coarseindex 
)

It computes the local contribution to K(c)

Definition at line 112 of file tpzdohrsubstruct.cpp.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_rc()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_rc ( TPZFMatrix< TVar > &  rc)

It computes the local contribution to r(c).

The method LoadWeightedResidual must be called before this one.

Definition at line 91 of file tpzdohrsubstruct.cpp.

References TPZMatrix< TVar >::Multiply().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_rc_local()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_rc_local ( TPZFMatrix< TVar > &  residual_local,
TPZFMatrix< TVar > &  rc_local 
) const

It computes the local contribution to r(c).

The method LoadWeightedResidual must be called before this one.

It computes the local contribution to r(c). The method LoadWeightedResidual must be called before this one.

Definition at line 105 of file tpzdohrsubstruct.cpp.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_v1()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_v1 ( TPZFMatrix< TVar > &  v1,
TPZFMatrix< TVar > &  invKc_rc 
)

It computes the local contribution to v1.

Parameters
v1
invKc_rcis the product K(c)_inverted*r(c) Of course r(c) must be computed, using Contribute_rc(), before calling this method

Definition at line 123 of file tpzdohrsubstruct.cpp.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_v1_local()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_v1_local ( TPZFMatrix< TVar > &  v1_local,
TPZFMatrix< TVar > &  invKc_rc 
) const

It computes the local contribution to v1.

Parameters
v1_local
invKc_rcis the product K(c)_inverted*r(c)Of course r(c) must be computed, using Contribute_rc(), before calling this method

Definition at line 149 of file tpzdohrsubstruct.cpp.

References TPZFMatrix< TVar >::Resize().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_v2()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_v2 ( TPZFMatrix< TVar > &  v2)

It computes the local contribution to v2.

Definition at line 156 of file tpzdohrsubstruct.cpp.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_v2_local()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_v2_local ( TPZFMatrix< TVar > &  residual_local,
TPZFMatrix< TVar > &  v2_local 
)

It computes the local contribution to v2.

Solving the system for zi

Definition at line 177 of file tpzdohrsubstruct.cpp.

References TPZMatrix< TVar >::Add(), TPZMatrix< TVar >::Cols(), LOGPZ_DEBUG, TPZMatrix< TVar >::Print(), TPZFMatrix< TVar >::Resize(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_v3()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_v3 ( TPZFMatrix< TVar > &  v3,
const TPZFMatrix< TVar > &  r,
TPZFMatrix< TVar > &  v1Plusv2 
) const

It computes the local contribution to v(3)

Parameters
v3
ris the global residual
v1Plusv2is the sum "v1 + v2"

Definition at line 255 of file tpzdohrsubstruct.cpp.

References TPZFMatrix< TVar >::GetVal(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Contribute_v3_local()

template<class TVar >
void TPZDohrSubstruct< TVar >::Contribute_v3_local ( TPZFMatrix< TVar > &  v3,
TPZFMatrix< TVar > &  v1Plusv2 
) const

It computes the local contribution to v(3)

Parameters
v3
v1Plusv2is the sum "v1 + v2"

Definition at line 298 of file tpzdohrsubstruct.cpp.

References TPZFMatrix< TVar >::Redim().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ContributeDiagonalLocal()

template<class TVar >
void TPZDohrSubstruct< TVar >::ContributeDiagonalLocal ( TPZFMatrix< TVar > &  StiffnessDiag)

Computes the contribution of each substructure node to global Stiffness diagonal (or something like that).

Parameters
StiffnessDiagis the diagonal of the stiffness matrix

Definition at line 493 of file tpzdohrsubstruct.cpp.

References LOGPZ_DEBUG, and TPZFMatrix< TVar >::Resize().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ContributeGlobalDiagonal()

template<class TVar >
void TPZDohrSubstruct< TVar >::ContributeGlobalDiagonal ( TPZFMatrix< TVar > &  StiffnessDiag)

Computes the contribution of each substructure node to global Stiffness diagonal (or something like that).

Parameters
StiffnessDiagis the diagonal of the stiffness matrix

Definition at line 454 of file tpzdohrsubstruct.cpp.

References LOGPZ_DEBUG.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ContributeKU()

template<class TVar >
void TPZDohrSubstruct< TVar >::ContributeKU ( const TVar  alpha,
const TPZFMatrix< TVar > &  uglobal,
TPZFMatrix< TVar > &  z 
) const

Contribute to the global matrix vector multiplication. It's needed to the MultAdd method of TPZDohrMatrix.

Definition at line 597 of file tpzdohrsubstruct.cpp.

References TPZMatrix< TVar >::Cols(), TPZMatrix< TVar >::Get(), LOGPZ_DEBUG, TPZMatrix< TVar >::Print(), TPZMatrix< TVar >::Rows(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ContributeKULocal()

template<class TVar >
void TPZDohrSubstruct< TVar >::ContributeKULocal ( const TVar  alpha,
const TPZFMatrix< TVar > &  u,
TPZFMatrix< TVar > &  z 
) const

Compute the multiplication of the local stiffness matrix with the vector u.

Definition at line 688 of file tpzdohrsubstruct.cpp.

References TPZMatrix< TVar >::Cols(), TPZMatrix< TVar >::Get(), LOGPZ_DEBUG, TPZMatrix< TVar >::Print(), TPZMatrix< TVar >::Rows(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ContributeResidual()

template<class TVar >
void TPZDohrSubstruct< TVar >::ContributeResidual ( TPZFMatrix< TVar > &  u,
TPZFMatrix< TVar > &  r 
)

ContributeResidual.

Parameters
uis the global solution
ris the residual of the solution "u" It makes $ u(i)=R(i)*u $, $ r(i)=F(i)-K(i)*u(i) $ and $ R(i)*r(novo)=R(i)*r+r(i) $ where r(novo) is the new global residual It puts r(novo) in r

Definition at line 61 of file tpzdohrsubstruct.cpp.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ContributeTestV1()

template<class TVar >
void TPZDohrSubstruct< TVar >::ContributeTestV1 ( TPZFMatrix< TVar > &  testV1,
int  NumCoarse 
)

◆ GetCoarseResidual()

template<class TVar>
void TPZDohrSubstruct< TVar >::GetCoarseResidual ( TPZFMatrix< TVar > &  rc,
TPZVec< int > &  indices 
)

Assembles the contribution to the coarse residual.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ GetCoarseStiffness()

template<class TVar>
void TPZDohrSubstruct< TVar >::GetCoarseStiffness ( TPZMatrix< TVar > &  stiff,
TPZVec< int > &  indices 
)

Assembles the coarse dof stiffness matrix.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Initialize()

template<class TVar >
void TPZDohrSubstruct< TVar >::Initialize ( )

Initializes the substructure.

Definition at line 758 of file tpzdohrsubstruct.cpp.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ LoadWeightedResidual()

template<class TVar >
void TPZDohrSubstruct< TVar >::LoadWeightedResidual ( const TPZFMatrix< TVar > &  r_global)

It computes W(i)*R(i)*r and stores it in fLocalWeightedResidual;.

Parameters
r_globalis the residual to the global solution

Definition at line 79 of file tpzdohrsubstruct.cpp.

References TPZFMatrix< TVar >::GetVal().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ NumEquations()

template<class TVar>
int TPZDohrSubstruct< TVar >::NumEquations ( ) const
inline

Definition at line 175 of file tpzdohrsubstruct.h.

References TPZDohrSubstruct< TVar >::fNEquations.

◆ PrepareSystems()

template<class TVar >
void TPZDohrSubstruct< TVar >::PrepareSystems ( )

It prepares the datas for solving systems for phi and zi.

Tou mantendo isso abaixo só pra fatoracao de K ser TVarizada e se obter os pivos nulos (caso haja)

Definition at line 766 of file tpzdohrsubstruct.cpp.

References ELU, LOGPZ_DEBUG, TPZFMatrix< TVar >::MultAdd(), TPZFMatrix< TVar >::Transpose(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ Print()

template<class TVar >
void TPZDohrSubstruct< TVar >::Print ( std::ostream &  out) const

Definition at line 318 of file tpzdohrsubstruct.cpp.

References EMathematicaInput.

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ ReallocMatRed()

template<class TVar>
void TPZDohrSubstruct< TVar >::ReallocMatRed ( )
inline

◆ SolveSystemPhi()

template<class TVar >
void TPZDohrSubstruct< TVar >::SolveSystemPhi ( )

Solves the system for Phi and for v2.

It stores the results in fPhiC and fzi

Acho q posso comentar essas duas linhas abaixo, pq naum tou vendo aplicacao pra temp2

Definition at line 354 of file tpzdohrsubstruct.cpp.

References LOGPZ_DEBUG, TPZMatrix< TVar >::Print(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

◆ SolveSystemZi()

template<class TVar >
void TPZDohrSubstruct< TVar >::SolveSystemZi ( )

Solves the system for zi.

It stores the results in fzi

Definition at line 398 of file tpzdohrsubstruct.cpp.

References TPZMatrix< TVar >::Add(), TPZFMatrix< TVar >::Resize(), and TPZFMatrix< TVar >::Zero().

Referenced by TPZDohrSubstruct< TVar >::ReallocMatRed().

Member Data Documentation

◆ fAdjustSolution

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fAdjustSolution

Solution vector which needs to be added to the converged system solution.

This variable is initialized and set in the AdjustResidual method

Definition at line 235 of file tpzdohrsubstruct.h.

◆ fBoundaryEqs

template<class TVar>
TPZVec<int> TPZDohrSubstruct< TVar >::fBoundaryEqs

Equations corresponding to boundary of the substructure.

Definition at line 224 of file tpzdohrsubstruct.h.

◆ fC

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fC

Constraint definition associated with this substructure.

Input data

Definition at line 186 of file tpzdohrsubstruct.h.

◆ fC_star

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fC_star

Variables needed to solve the systems for phi and zi.

Definition at line 167 of file tpzdohrsubstruct.h.

◆ fCoarseIndex

template<class TVar>
TPZVec<int> TPZDohrSubstruct< TVar >::fCoarseIndex

Global index associated with each coarse degree of freedom.

Input data

Definition at line 194 of file tpzdohrsubstruct.h.

◆ fCoarseNodes

template<class TVar>
TPZVec<int> TPZDohrSubstruct< TVar >::fCoarseNodes

Definition at line 182 of file tpzdohrsubstruct.h.

◆ fGlobalEqs

template<class TVar>
TPZVec<std::pair<int,int> > TPZDohrSubstruct< TVar >::fGlobalEqs

Indices of the equations in the global system corresponding to the local equations.

Definition at line 210 of file tpzdohrsubstruct.h.

◆ fGlobalIndex

template<class TVar>
TPZVec<int> TPZDohrSubstruct< TVar >::fGlobalIndex

Global vector indices of the equations/degrees of freedom.

Input data

Definition at line 206 of file tpzdohrsubstruct.h.

◆ fInternalEqs

template<class TVar>
TPZVec<int> TPZDohrSubstruct< TVar >::fInternalEqs

Internal nodes.

Data

Definition at line 222 of file tpzdohrsubstruct.h.

◆ finv

template<class TVar>
TPZStepSolver<TVar> TPZDohrSubstruct< TVar >::finv

Definition at line 169 of file tpzdohrsubstruct.h.

◆ fInvertedInternalStiffness

template<class TVar>
TPZStepSolver<TVar> TPZDohrSubstruct< TVar >::fInvertedInternalStiffness
mutable

Inverted (LU or Cholesky or LDLt) stiffness matrix for the internal degrees of freedom.

Computed

Definition at line 216 of file tpzdohrsubstruct.h.

◆ fInvertedStiffness

template<class TVar>
TPZStepSolver<TVar> TPZDohrSubstruct< TVar >::fInvertedStiffness

Inverted (LU or Cholesky) stiffness matrix.

Computed

Definition at line 213 of file tpzdohrsubstruct.h.

◆ fKCi

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fKCi

Stiffness matrix associated with the constraints.

Computed

Definition at line 203 of file tpzdohrsubstruct.h.

◆ fKCInvert

template<class TVar>
TPZStepSolver<TVar> TPZDohrSubstruct< TVar >::fKCInvert

Inverted restraint matrix.

Computed

Definition at line 219 of file tpzdohrsubstruct.h.

◆ fKeC_star

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fKeC_star

Definition at line 168 of file tpzdohrsubstruct.h.

◆ fLocalLoad

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fLocalLoad

Local load vector.

Data

Definition at line 227 of file tpzdohrsubstruct.h.

◆ fLocalWeightedResidual

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fLocalWeightedResidual

Local weighted residual - $ W(i)*R(i)*r $.

Definition at line 229 of file tpzdohrsubstruct.h.

◆ fNEquations

template<class TVar>
int TPZDohrSubstruct< TVar >::fNEquations

Number of equations of the substructure.

Definition at line 173 of file tpzdohrsubstruct.h.

Referenced by TPZDohrSubstruct< TVar >::NumEquations().

◆ fNullPivots

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fNullPivots

Definition at line 170 of file tpzdohrsubstruct.h.

◆ fPhiC

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fPhiC

Vectors associated with each constraint.

Definition at line 188 of file tpzdohrsubstruct.h.

◆ fPhiC_Weighted_Condensed

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fPhiC_Weighted_Condensed

Phi * W matrix and condensed to the equations which are part of the global system.

Definition at line 190 of file tpzdohrsubstruct.h.

◆ fStiffness

template<class TVar>
TPZAutoPointer<TPZMatrix<TVar> > TPZDohrSubstruct< TVar >::fStiffness

Stiffness matrix associated with the substructure.

Input data

Definition at line 200 of file tpzdohrsubstruct.h.

◆ fWeights

template<class TVar>
TPZManVector<TVar> TPZDohrSubstruct< TVar >::fWeights

Weights associated with each variable/equation.

Computed

Definition at line 197 of file tpzdohrsubstruct.h.

◆ fWeightType

template<class TVar>
TPZDohrSubstruct< TVar >::EWeightType TPZDohrSubstruct< TVar >::fWeightType
static

Definition at line 41 of file tpzdohrsubstruct.h.

◆ fzi

template<class TVar>
TPZFMatrix<TVar> TPZDohrSubstruct< TVar >::fzi

Needed to compute v2. Is the solution of a system.

Definition at line 231 of file tpzdohrsubstruct.h.


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