OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1_2d_supg.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 - 2013 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 tr1_2d_supg_h
36 #define tr1_2d_supg_h
37 
38 #include "supgelement.h"
39 #include "primaryfield.h"
40 #include "spatiallocalizer.h"
41 #include "zznodalrecoverymodel.h"
43 #include "sprnodalrecoverymodel.h"
44 #include "leplic.h"
45 #include "levelsetpcs.h"
46 
48 
49 #define _IFT_TR1_2D_SUPG_Name "tr1supg"
50 #define _IFT_Tr1SUPG_pvof "pvof"
51 #define _IFT_Tr1SUPG_vof "vof"
52 #define _IFT_Tr1SUPG2_mat0 "mat0"
53 #define _IFT_Tr1SUPG2_mat1 "mat1"
54 
55 
56 namespace oofem {
57 class FEI2dTrLin;
58 
68 class TR1_2D_SUPG : public SUPGElement,
72 {
73 protected:
75 
76  //double a[3];
77  double b [ 3 ];
78  double c [ 3 ];
79  double area;
80 
81 public:
82  TR1_2D_SUPG(int n, Domain * d);
83  virtual ~TR1_2D_SUPG();
84 
85  virtual FEInterpolation *giveInterpolation() const;
86 
87  virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep);
88  virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep);
89  virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep);
90  virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep);
91  virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep);
92  virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep);
93  virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep);
94  virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep);
95  virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep);
96  virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep);
97  virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep) {
98  answer.resize(3, 6);
99  answer.zero();
100  }
101  virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep) {
102  answer.resize(3);
103  answer.zero();
104  }
105  virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep);
106  virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep);
107  virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep);
108  virtual void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep);
109  virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep);
110 
111  virtual void computeSlipWithFrictionBCTerm_MB(FloatMatrix &answer, Load *load, int side, TimeStep *tStep);
112  virtual void computePenetrationWithResistanceBCTerm_MB(FloatMatrix &answer, Load *load, int side, TimeStep *tStep);
113  virtual void computeOutFlowBCTerm_MB(FloatMatrix &answer, int side, TimeStep *tStep);
114 
115  virtual void computeHomogenizedReinforceTerm_MB(FloatMatrix &answer, Load *load, TimeStep *tStep);
116  virtual void computeHomogenizedReinforceTerm_MC(FloatMatrix &answer, Load *load, TimeStep *tStep);
117 
118  virtual void updateStabilizationCoeffs(TimeStep *tStep);
119  virtual double computeCriticalTimeStep(TimeStep *tStep);
120 
121  // definition
122  virtual const char *giveClassName() const { return "TR1_2D_SUPG"; }
123  virtual const char *giveInputRecordName() const { return _IFT_TR1_2D_SUPG_Name; }
124  virtual MaterialMode giveMaterialMode() { return _2dFlow; }
125 
126  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
127  virtual int computeNumberOfDofs();
129  virtual void giveInputRecord(DynamicInputRecord &input);
130  virtual void updateYourself(TimeStep *tStep);
132  virtual int checkConsistency();
133 
134  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
135  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
136 
138 
140  const FloatArray &coords, IntArray &dofId, ValueModeType mode,
141  TimeStep *tStep);
142 
143  virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag);
144  virtual void formMaterialVolumePoly(Polygon &matvolpoly, LEPlic *matInterface,
145  const FloatArray &normal, const double p, bool updFlag);
146  virtual void formVolumeInterfacePoly(Polygon &matvolpoly, LEPlic *matInterface,
147  const FloatArray &normal, const double p, bool updFlag);
148  virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume);
149  virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool updFlag);
150  virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag);
151  virtual Element *giveElement() { return this; }
152  virtual double computeMyVolume(LEPlic *matInterface, bool updFlag);
153  virtual double computeVolumeAround(GaussPoint *gp);
154  virtual double computeCriticalLEPlicTimeStep(TimeStep *tStep);
155 
157  InternalStateType type, TimeStep *tStep);
158 
160  virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap);
161  virtual int SPRNodalRecoveryMI_giveNumberOfIP();
163 
164  virtual double LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep);
165  virtual void LS_PCS_computedN(FloatMatrix &answer);
166  virtual double LS_PCS_computeVolume() { return area; }
167  virtual double LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep);
168  virtual void LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi);
169 
170  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
171 
172 #ifdef __OOFEG
174  int node, TimeStep *tStep);
175  // Graphics output
176  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
177  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
178  //virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) {}
179 #endif
180 
181  virtual void printOutputAt(FILE *file, TimeStep *tStep);
182 
183 protected:
184  virtual void giveLocalVelocityDofMap(IntArray &map);
185  virtual void giveLocalPressureDofMap(IntArray &map);
186 
187  void computeNMtrx(FloatArray &answer, GaussPoint *gp);
188  virtual void computeGaussPoints();
189 
190  virtual void computeDeviatoricStrain(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
191  virtual void initGeometry();
192 };
193 } // end namespace oofem
194 #endif // tr1_2d_supg_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms (generalized mass matrix with stabilization terms) for momentum balance e...
Definition: tr1_2d_supg.C:147
void computeNMtrx(FloatArray &answer, GaussPoint *gp)
Definition: tr1_2d_supg.C:1544
virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for momentum balance equations(s).
Definition: tr1_2d_supg.C:323
virtual int EIPrimaryFieldI_evaluateFieldVectorAt(FloatArray &answer, PrimaryField &pf, const FloatArray &coords, IntArray &dofId, ValueModeType mode, TimeStep *tStep)
Evaluates the value of field at given point of interest (should be located inside receiver's volume) ...
Definition: tr1_2d_supg.C:1867
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: tr1_2d_supg.C:1808
The element interface required by ZZNodalRecoveryModel.
virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes nonlinear advection terms for momentum balance equations(s).
Definition: tr1_2d_supg.C:205
Element interface for LEPlic class representing Lagrangian-Eulerian (moving) material interface...
Definition: leplic.h:58
Class and object Domain.
Definition: domain.h:115
virtual void computeHomogenizedReinforceTerm_MC(FloatMatrix &answer, Load *load, TimeStep *tStep)
Definition: tr1_2d_supg.C:846
#define _IFT_TR1_2D_SUPG_Name
Definition: tr1_2d_supg.h:49
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr1_2d_supg.C:135
virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes advection terms for mass conservation equation.
Definition: tr1_2d_supg.C:473
virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the linear advection term for mass conservation equation.
Definition: tr1_2d_supg.C:456
The element interface required by ZZNodalRecoveryModel.
Abstract class representing field of primary variables (those, which are unknown and are typically as...
Definition: primaryfield.h:104
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
Definition: tr1_2d_supg.C:1011
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
Class implementing element body load, acting over whole element volume (e.g., the dead weight)...
Definition: bodyload.h:49
Element interface for LevelSetPCS class representing level-set like material interface.
Definition: levelsetpcs.h:68
virtual int checkConsistency()
Used to check consistency and initialize some element geometry data (area,b,c).
Definition: tr1_2d_supg.C:1538
virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume)
Truncates given material polygon to receiver.
Definition: tr1_2d_supg.C:1749
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr1_2d_supg.C:98
virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag)
Assembles receiver volume.
Definition: tr1_2d_supg.C:1767
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
Abstract base class for all finite elements.
Definition: element.h:145
virtual void giveLocalPressureDofMap(IntArray &map)
Definition: tr1_2d_supg.C:2260
int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
Definition: tr1_2d_supg.C:2268
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr1_2d_supg.C:2317
General stabilized SUPG/PSPG element for CFD analysis.
Definition: supgelement.h:58
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes diffusion derivative terms for mass conservation equation.
Definition: tr1_2d_supg.h:97
Class implementing an array of integers.
Definition: intarray.h:61
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr1_2d_supg.C:1938
MatResponseMode
Describes the character of characteristic material matrix.
virtual void computeHomogenizedReinforceTerm_MB(FloatMatrix &answer, Load *load, TimeStep *tStep)
Definition: tr1_2d_supg.C:818
virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for mass conservation equation.
Definition: tr1_2d_supg.C:557
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Definition: tr1_2d_supg.C:1466
virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool updFlag)
Computes the receiver center (in updated Lagrangian configuration).
Definition: tr1_2d_supg.C:1845
virtual double LS_PCS_computeVolume()
Returns receiver's volume.
Definition: tr1_2d_supg.h:166
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Element interface class.
Definition: primaryfield.h:58
virtual const char * giveClassName() const
Definition: tr1_2d_supg.h:122
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: tr1_2d_supg.C:1917
Class representing a 2d triangular linear interpolation based on area coordinates.
Definition: fei2dtrlin.h:44
virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
Definition: tr1_2d_supg.C:867
virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
Computes corresponding volume fraction to given interface position.
Definition: tr1_2d_supg.C:1594
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Updates the stabilization coefficients used for CBS and SUPG algorithms.
Definition: tr1_2d_supg.C:1089
virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms for mass conservation equation.
Definition: tr1_2d_supg.C:541
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr1_2d_supg.C:1910
virtual double LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep)
Evaluates F in level set equation of the form where for interface position driven by flow with speed...
Definition: tr1_2d_supg.C:2033
virtual void initGeometry()
Definition: tr1_2d_supg.C:1504
virtual void computeSlipWithFrictionBCTerm_MB(FloatMatrix &answer, Load *load, int side, TimeStep *tStep)
Computes Lhs term due to applied slip with friction bc.
Definition: tr1_2d_supg.C:573
virtual void giveLocalVelocityDofMap(IntArray &map)
Definition: tr1_2d_supg.C:2254
virtual void formMaterialVolumePoly(Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag)
Assembles the true element material polygon (takes receiver vof into accout).
Definition: tr1_2d_supg.C:1608
virtual void LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi)
Returns VOF fractions for each material on element according to nodal values of level set function (p...
Definition: tr1_2d_supg.C:2165
virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for mass conservation equation.
Definition: tr1_2d_supg.h:101
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: tr1_2d_supg.C:1990
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: tr1_2d_supg.C:2011
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: tr1_2d_supg.C:1955
virtual void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
Definition: tr1_2d_supg.C:966
virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for mass conservation equation with respect to nodal veloc...
Definition: tr1_2d_supg.C:500
virtual double computeCriticalTimeStep(TimeStep *tStep)
Computes the critical time increment.
Definition: tr1_2d_supg.C:1459
Abstract base class representing Lagrangian-Eulerian (moving) material interfaces.
Definition: leplic.h:153
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr1_2d_supg.C:2287
virtual void LS_PCS_computedN(FloatMatrix &answer)
Returns gradient of shape functions.
Definition: tr1_2d_supg.C:2153
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void computeDeviatoricStrain(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Definition: tr1_2d_supg.C:1489
Class representing 2D polygon.
Definition: geotoolbox.h:91
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
TR1_2D_SUPG(int n, Domain *d)
Definition: tr1_2d_supg.C:73
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: tr1_2d_supg.h:124
virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for momentum balance equations(s) with respect to nodal ve...
Definition: tr1_2d_supg.C:249
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
Abstract base class representing Level Set representation of material interfaces. ...
Definition: levelsetpcs.h:114
static FEI2dTrLin interp
Definition: tr1_2d_supg.h:74
CharType
Definition: chartype.h:87
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: tr1_2d_supg.C:83
Class Interface.
Definition: interface.h:82
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
Class representing the a dynamic Input Record.
virtual double LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep)
Evaluates S in level set equation of the form where .
Definition: tr1_2d_supg.C:2054
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Definition: tr1_2d_supg.C:1968
The spatial localizer element interface associated to spatial localizer.
virtual FEInterpolation * giveInterpolation() const
Definition: tr1_2d_supg.C:95
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
virtual double computeMyVolume(LEPlic *matInterface, bool updFlag)
Computes the volume of receiver.
Definition: tr1_2d_supg.C:1790
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
Definition: tr1_2d_supg.C:1972
virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Computes the derivative of diffusion terms for momentum balance equations(s) with respect to nodal ve...
Definition: tr1_2d_supg.C:351
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual Element * giveElement()
Return number of receiver's element.
Definition: tr1_2d_supg.h:151
Class representing 2d linear triangular element for solving incompressible fluid with SUPG solver...
Definition: tr1_2d_supg.h:68
virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes SLIC stabilization term for momentum balance equation(s).
Definition: tr1_2d_supg.C:437
virtual const char * giveInputRecordName() const
Definition: tr1_2d_supg.h:123
Load is base abstract class for all loads.
Definition: load.h:61
double p
Line constant of line segment representing interface.
Definition: leplic.h:65
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual ~TR1_2D_SUPG()
Definition: tr1_2d_supg.C:79
virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for momentum balance equations(s).
Definition: tr1_2d_supg.C:385
Class representing integration point in finite element program.
Definition: gausspoint.h:93
FloatArray normal
Interface segment normal.
Definition: leplic.h:67
Class representing solution step.
Definition: timestep.h:80
virtual void formVolumeInterfacePoly(Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag)
Assembles receiver material polygon based solely on given interface line.
Definition: tr1_2d_supg.C:1640
virtual void computePenetrationWithResistanceBCTerm_MB(FloatMatrix &answer, Load *load, int side, TimeStep *tStep)
Computes Lhs contribution due to applied Penetration bc.
Definition: tr1_2d_supg.C:667
InternalStateMode
Determines the mode of internal variable.
virtual void computeOutFlowBCTerm_MB(FloatMatrix &answer, int side, TimeStep *tStep)
Computes Lhs contribution due to outflow BC.
Definition: tr1_2d_supg.C:759
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: tr1_2d_supg.C:1946
virtual double computeCriticalLEPlicTimeStep(TimeStep *tStep)
Computes critical time step.
Definition: tr1_2d_supg.C:1818
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: tr1_2d_supg.C:1980
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr1_2d_supg.C:89
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: tr1_2d_supg.C:123

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:32 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011