61 IntArray
Tr21Stokes :: momentum_ordering = {1, 2, 4, 5, 7, 8, 10, 11, 12, 13, 14, 15};
99 answer = {V_u, V_v, P_f};
115 if ( mtrx == ExternalForcesVector ) {
117 }
else if ( mtrx == InternalForcesVector ) {
120 OOFEM_ERROR(
"Unknown Type of characteristic mtrx.");
128 if ( mtrx == TangentStiffnessMatrix ) {
131 OOFEM_ERROR(
"Unknown Type of characteristic mtrx.");
138 FloatArray a_pressure, a_velocity, devStress, epsp, Nh, dNv(12);
139 double r_vol, pressure;
149 const FloatArray &lcoords = gp->giveNaturalCoordinates();
153 double dA = detJ * gp->giveWeight();
156 dNv(k) = B(0, k) = B(2, k + 1) = dN(j, 0);
157 dNv(k + 1) = B(1, k + 1) = B(2, k) = dN(j, 1);
166 momentum.
add(-pressure * dA, dNv);
167 conservation.
add(r_vol * dA, Nh);
183 for (
int i = 1; i <= nLoads; i++ ) {
197 for (
int i = 1; i <= nLoads; i++ ) {
199 if ((bload = dynamic_cast<BodyLoad*>(load))) {
211 if ( type != ExternalForcesVector ) {
223 const FloatArray &lcoords = gp->giveNaturalCoordinates();
225 double rho = mat->
give(
'd', gp);
227 double dA = detJ * gp->giveWeight();
230 for (
int j = 0; j < 6; j++ ) {
231 temparray(2 * j) +=
N(j) * rho * gVector(0) * dA;
232 temparray(2 * j + 1) +=
N(j) * rho * gVector(1) * dA;
244 if ( type != ExternalForcesVector ) {
251 int numberOfEdgeIPs = ( int ) ceil( ( load->
giveApproxOrder() + 1. ) / 2. ) * 2;
260 const FloatArray &lcoords = gp->giveNaturalCoordinates();
264 double dS = gp->giveWeight() * detJ;
275 for (
int j = 0; j < 3; j++ ) {
276 f(2 * j) +=
N(j) * t(0) * dS;
277 f(2 * j + 1) +=
N(j) * t(1) * dS;
293 FloatMatrix B(3, 12), EdB, K, G, Dp, DvT, C, Ed, dN;
294 FloatArray dN_V(12), Nlin, Ep, Cd, tmpA, tmpB;
301 const FloatArray &lcoords = gp->giveNaturalCoordinates();
305 double dA = detJ * gp->giveWeight();
307 for (
int j = 0, k = 0; j < 6; j++, k += 2 ) {
308 dN_V(k) = B(0, k) = B(2, k + 1) = dN(j, 0);
309 dN_V(k + 1) = B(1, k + 1) = B(2, k) = dN(j, 1);
392 for (
int i = 1; i <= n.
giveSize(); i++ ) {
393 answer(0) += n.
at(i) * velocities.
at(i*2-1);
394 answer(1) += n.
at(i) * velocities.
at(i*2);
397 for (
int i = 1; i <= n_lin.
giveSize(); i++ ) {
398 answer(2) += n_lin.
at(i) * pressures.
at(i);
405 if ( type == IST_Pressure ) {
407 if ( node == 1 || node == 2 || node == 3 ) {
414 }
else if ( node == 5 ) {
422 answer.
at(1) = ( a + b ) / 2;
CrossSection * giveCrossSection()
Tr21Stokes(int n, Domain *d)
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
problemScale giveProblemScale()
Returns scale in multiscale simulation.
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Abstract base class for all fluid materials.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Gives the dof ID mask for the element.
static FEI2dTrQuad interpolation_quad
Interpolation for geometry and velocity.
The element interface required by ZZNodalRecoveryModel.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
bcGeomType
Type representing the geometric character of loading.
double & at(int i)
Coefficient access function.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Class implementing element body load, acting over whole element volume (e.g., the dead weight)...
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep)
void clear()
Clears receiver (zero size).
virtual void computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
Computes the contribution of the given load at the given boundary surface in global coordinate system...
virtual double giveUnknown(ValueModeType mode, TimeStep *tStep)=0
The key method of class Dof.
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
Adds to the receiver the product .
virtual void giveStiffnessMatrices(FloatMatrix &dsdd, FloatArray &dsdp, FloatArray &dedd, double &dedp, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the 4 tangents of the compressible material response in 3D.
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
MatResponseMode
Describes the character of characteristic material matrix.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Class representing a general abstraction for finite element interpolation class.
virtual void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property 'aProperty'.
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
static IntArray momentum_ordering
Ordering of dofs in element. Used to assemble the element stiffness.
virtual void computeDeviatoricStressVector(FloatArray &stress_dev, double &epsp_vol, GaussPoint *gp, const FloatArray &eps, double pressure, TimeStep *tStep)
Computes the deviatoric stress vector and volumetric strain rate from given deviatoric strain and pre...
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
virtual void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
DofIDItem
Type representing particular dof type.
virtual double giveWeight()
Returns integration weight of receiver.
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Neumann type (prescribed flux).
static FEI2dTrLin interpolation_lin
Interpolation for pressure.
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation Jacobian on the requested boundary.
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Wrapper around element definition to provide FEICellGeometry interface.
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
IntArray bodyLoadArray
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver.
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
virtual int giveApproxOrder()=0
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
virtual bcType giveType() const
Returns receiver load type.
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
int numberOfGaussPoints
Number of integration points as specified by nip.
Class representing vector of real numbers.
void plusDyadSymmUpper(const FloatArray &a, double dV)
Adds to the receiver the dyadic product .
This abstract class represent a general base element class for fluid dynamic problems.
Implementation of matrix containing floating point numbers.
void assemble(const FloatArray &fe, const IntArray &loc)
Assembles the array fe (typically, the load vector of a finite element) into the receiver, using loc as location array.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Dof * giveDofWithID(int dofID) const
Returns DOF with given dofID; issues error if not present.
void zero()
Zeroes all coefficients of receiver.
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual Interface * giveInterface(InterfaceType t)
Interface requesting service.
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
The spatial localizer element interface associated to spatial localizer.
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
void zero()
Zeroes all coefficient of receiver.
InterfaceType
Enumerative type, used to identify interface type.
void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
Load is base abstract class for all loads.
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
virtual FEInterpolation * giveInterpolation() const
int giveSize() const
Returns the size of receiver.
the oofem namespace is to define a context or scope in which all oofem names are defined.
void assemble(const FloatMatrix &src, const IntArray &loc)
Assembles the contribution using localization array into receiver.
virtual bcValType giveBCValType() const
Returns receiver load type.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Maps the local boundary coordinates to global.
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)
Computes components values of load at given point - global coordinates (coordinates given)...
int giveNumberOfRows() const
Returns number of rows of receiver.
Load * giveLoad(int n)
Service for accessing particular domain load.
virtual int SetUpPointsOnLine(int nPoints, MaterialMode mode)
Sets up receiver's integration points on unit line integration domain.
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
virtual void edgeEvalN(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of edge interpolation functions (shape functions) at given point.
Class representing integration point in finite element program.
IntArray boundaryLoadArray
Class representing solution step.
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
int numberOfDofMans
Number of dofmanagers.
void add(const FloatArray &src)
Adds array src to receiver.
virtual void computeField(ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Class representing Gaussian-quadrature integration rule.
static IntArray edge_ordering[3]
Ordering of dofs on edges. Used to assemble edge loads.
void resize(int s)
Resizes receiver towards requested size.
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
static IntArray conservation_ordering