OOFEM 3.0
Loading...
Searching...
No Matches
structuralelementevaluator.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#ifndef structuralelementevaluator_h
36#define structuralelementevaluator_h
37
38#include "iga/iga.h"
39#include "matresponsemode.h"
40
41namespace oofem {
57{
58protected:
62
65 virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep);
66 virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep);
67
72 virtual IntegrationRule *giveMassMtrxIntegrationRule() { return nullptr; }
81 virtual void giveMassMtrxIntegrationMask(IntArray &answer) { answer.clear(); }
92 virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
105 virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass);
106
107protected:
108 virtual Element *giveElement() = 0;
109
113 virtual void computeNMatrixAt(FloatMatrix &answer, GaussPoint *gp) = 0;
114 virtual void computeBMatrixAt(FloatMatrix &answer, GaussPoint *gp) = 0;
115 virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
116 virtual double computeVolumeAround(GaussPoint *gp) { return 0.; }
117 virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, bool useUpdatedGpRecord = false);
118 void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer) {
119 this->giveElement()->computeVectorOf(u, tStep, answer);
120 }
121 void computeVectorOf(PrimaryField &field, const IntArray &dofIdMask, ValueModeType u, TimeStep *tStep, FloatArray &answer) {
122 this->giveElement()->computeVectorOf(field, dofIdMask, u, tStep, answer);
123 }
124 bool isActivated(TimeStep *tStep) { return true; }
125 void updateInternalState(TimeStep *tStep);
126
134 virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) = 0;
142 virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) = 0;
150 // virtual void computeConstitutiveMatrix_dPdF_At(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) = 0;
151
159 void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, FloatArray &u);
160
161 /*
162 * Assembles the code numbers of given integration element (sub-patch)
163 * This is done by obtaining list of nonzero shape functions and
164 * by collecting the code numbers of nodes corresponding to these
165 * shape functions
166 * @returns returns nonzero if integration rule code numbers differ from element code numbers
167 */
168 //virtual int giveIntegrationElementCodeNumbers(IntArray &answer, Element *elem,
169 // IntegrationRule *ie,);
170
179 IntegrationRule *ie);
180#ifdef __OOFEG
182#endif
183};
184} // end namespace oofem
185#endif //structuralelementevaluator_h
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, bool useUpdatedGpRecord=false)
virtual void computeNMatrixAt(FloatMatrix &answer, GaussPoint *gp)=0
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep)
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass)
virtual double computeVolumeAround(GaussPoint *gp)
int rotationMatrixDefined
Flag indicating if transformation matrix has been already computed.
virtual void giveMassMtrxIntegrationMask(IntArray &answer)
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, FloatArray &u)
virtual Element * giveElement()=0
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
void computeVectorOf(PrimaryField &field, const IntArray &dofIdMask, ValueModeType u, TimeStep *tStep, FloatArray &answer)
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
friend void drawIGAPatchDeformedGeometry(Element *elem, StructuralElementEvaluator *se, oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition iga.C:1001
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
virtual void computeBMatrixAt(FloatMatrix &answer, GaussPoint *gp)=0
virtual int giveIntegrationElementLocalCodeNumbers(IntArray &answer, Element *elem, IntegrationRule *ie)
virtual IntegrationRule * giveMassMtrxIntegrationRule()
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]

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