OOFEM 3.0
Loading...
Searching...
No Matches
pfemelement.h
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35// ********************************************************************************
36// *** GENERAL ELEMENT CLASS FOR FLUID DYNAMIC PROBLEMS SOLVED WITH PFEM METHOD ***
37// ********************************************************************************
38
39#ifndef pfemelement_h
40#define pfemelement_h
41
42
43#include "fmelement.h"
44#include "femcmpnn.h"
45#include "domain.h"
46#include "floatmatrix.h"
47#include "material.h"
48
49#include "primaryfield.h"
50
51namespace oofem {
52class TimeStep;
53class Node;
54class Material;
55class GaussPoint;
56class FloatMatrix;
57class FloatArray;
58class IntArray;
59
66class PFEMElement : public FMElement
67{
68public:
70 PFEMElement(int, Domain *);
72 virtual ~PFEMElement();
73
75 void initializeFrom(InputRecord &ir, int priority) override;
76
77 // characteristic matrix
78 void giveCharacteristicMatrix(FloatMatrix & answer, CharType, TimeStep *) override;
79 void giveCharacteristicVector(FloatArray & answer, CharType, ValueModeType, TimeStep *) override;
80
82 virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *) = 0;
84 virtual void computeDiagonalMassMtrx(FloatMatrix &answer, TimeStep *) = 0;
86 virtual double computeCriticalTimeStep(TimeStep *tStep) = 0;
87
88 void updateInternalState(TimeStep *tStep) override;
89 void printOutputAt(FILE *file, TimeStep *tStep) override;
90 int checkConsistency() override;
91
96
97 void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) override;
98
99 // definition
100 const char *giveClassName() const override { return "PFEMElement"; }
101
103 virtual const IntArray &giveVelocityDofMask() const = 0;
105 virtual const IntArray &givePressureDofMask() const = 0;
106
107#ifdef __OOFEG
109 int node, TimeStep *atTime);
110 //
111 // Graphics output
112 //
113 //virtual void drawYourself(oofegGraphicContext&);
114 //virtual void drawRawGeometry(oofegGraphicContext&) {}
115 //virtual void drawDeformedGeometry(oofegGraphicContext&, UnknownType) {}
116#endif
117
127 virtual void computeBodyLoadVectorAt(FloatArray &answer, BodyLoad *load, TimeStep *tStep, ValueModeType mode) = 0;
128
129
131 virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *) = 0;
133 virtual void computeDeviatoricStressDivergence(FloatArray &answer, TimeStep *atTime) = 0;
135 virtual void computeBMatrix(FloatMatrix &answer, GaussPoint *gp) = 0;
137 virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *atTime) = 0; //K
139 virtual void computePressureLaplacianMatrix(FloatMatrix &answer, TimeStep *atTime) = 0;
141 virtual void computeGradientMatrix(FloatMatrix &answer, TimeStep *atTime) = 0; //G
143 virtual void computeDivergenceMatrix(FloatMatrix &answer, TimeStep *atTime) = 0; //D
145 virtual void computePrescribedRhsVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode) = 0;
146};
147} // end namespace oofem
148#endif // pfemelement_h
FMElement(int n, Domain *aDomain)
Definition fmelement.C:38
virtual void computeDiagonalMassMtrx(FloatMatrix &answer, TimeStep *)=0
virtual void computeDeviatoricStressDivergence(FloatArray &answer, TimeStep *atTime)=0
Calculates the divergence of the deviatoric stress.
virtual FEInterpolation * giveVelocityInterpolation()=0
Returns the interpolation for velocity.
virtual double computeCriticalTimeStep(TimeStep *tStep)=0
Calculates critical time step.
virtual void computePrescribedRhsVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)=0
Calculates the prescribed velocity vector for the right hand side of the pressure equation.
void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *) override
Definition pfemelement.C:76
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *atTime)=0
Calculates the stiffness matrix.
virtual const IntArray & givePressureDofMask() const =0
Returns mask of pressure Dofs.
virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *)=0
virtual ~PFEMElement()
Destructor.
Definition pfemelement.C:63
virtual void computeDivergenceMatrix(FloatMatrix &answer, TimeStep *atTime)=0
Calculates the velocity divergence matrix.
virtual FEInterpolation * givePressureInterpolation()=0
Returns the interpolation for the pressure.
void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) override
virtual void computePressureLaplacianMatrix(FloatMatrix &answer, TimeStep *atTime)=0
Calculates the pressure laplacian matrix.
virtual const IntArray & giveVelocityDofMask() const =0
Returns mask of velocity Dofs.
PFEMElement(int, Domain *)
Constructor.
Definition pfemelement.C:58
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *)=0
Computes deviatoric stress vector in given integration point and solution step from given total strai...
void updateInternalState(TimeStep *tStep) override
void printOutputAt(FILE *file, TimeStep *tStep) override
virtual void computeGradientMatrix(FloatMatrix &answer, TimeStep *atTime)=0
Calculates the pressure gradient matrix.
void giveCharacteristicVector(FloatArray &answer, CharType, ValueModeType, TimeStep *) override
virtual void computeBodyLoadVectorAt(FloatArray &answer, BodyLoad *load, TimeStep *tStep, ValueModeType mode)=0
virtual void computeBMatrix(FloatMatrix &answer, GaussPoint *gp)=0
Calculates the element shape function derivative matrix.
virtual int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *atTime)
int checkConsistency() override
const char * giveClassName() const override
void initializeFrom(InputRecord &ir, int priority) override
Initializes receiver acording to object description stored in input record.
Definition pfemelement.C:67
InternalStateMode
Determines the mode of internal variable.

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011