96 answer = {V_u, V_v, P_f};
112 return this->
bubble.get();
125 if ( mtrx == ExternalForcesVector ) {
127 }
else if ( mtrx == InternalForcesVector ) {
130 OOFEM_ERROR(
"Unknown Type of characteristic mtrx.");
137 if ( mtrx == TangentStiffnessMatrix ) {
140 OOFEM_ERROR(
"Unknown Type of characteristic mtrx.");
147 FloatArray a_pressure, a_velocity, devStress, epsp,
N, dNv(8);
148 double r_vol, pressure;
158 const FloatArray &lcoords = gp->giveNaturalCoordinates();
162 double dA = detJ * gp->giveWeight();
164 for (
int j = 0, k = 0; j < 3; j++, k += 2 ) {
165 dNv(k) = B(0, k) = B(2, k + 1) = dN(j, 0);
166 dNv(k + 1) = B(1, k + 1) = B(2, k) = dN(j, 1);
170 dNv(6) = B(0, 6) = B(2, 7) = 27. * ( dN(0, 0) *
N(1) *
N(2) +
N(0) * dN(1, 0) *
N(2) +
N(0) *
N(1) * dN(2, 0) );
171 dNv(7) = B(1, 7) = B(2, 6) = 27. * ( dN(0, 1) *
N(1) *
N(2) +
N(0) * dN(1, 1) *
N(2) +
N(0) *
N(1) * dN(2, 1) );
179 momentum.
add(-pressure * dA, dNv);
180 conservation.
add(r_vol * dA, N);
191 int load_number, load_id;
200 for (
int i = 1; i <= nLoads; i++ ) {
213 for (
int i = 1; i <= nLoads; i++ ) {
215 if ((bLoad = dynamic_cast<BodyLoad*>(load))) {
228 if ( type != ExternalForcesVector ) {
239 const FloatArray &lcoords = gp->giveNaturalCoordinates();
243 double dA = detJ * gp->giveWeight();
246 for (
int j = 0; j < 3; j++ ) {
247 temparray(2 * j) +=
N(j) * rho * gVector(0) * dA;
248 temparray(2 * j + 1) +=
N(j) * rho * gVector(1) * dA;
251 temparray(6) +=
N(0) *
N(1) *
N(2) * rho * gVector(0) * dA;
252 temparray(7) +=
N(0) *
N(1) *
N(2) * rho * gVector(1) * dA;
263 if ( type != ExternalForcesVector ) {
270 int numberOfEdgeIPs = ( int ) ceil( ( load->
giveApproxOrder() + 2. ) / 2. );
279 const FloatArray &lcoords = gp->giveNaturalCoordinates();
283 double dS = gp->giveWeight() * detJ;
294 for (
int j = 0; j < 2; j++ ) {
295 f(2 * j) +=
N(j) * t(0) * dS;
296 f(2 * j + 1) +=
N(j) * t(1) * dS;
312 FloatMatrix B(3, 8), EdB, K, G, Dp, DvT, C, Ed, dN;
319 const FloatArray &lcoords = gp->giveNaturalCoordinates();
322 double dA = detJ * gp->giveWeight();
324 for (
int j = 0, k = 0; j < 3; j++, k += 2 ) {
325 dNv(k) = B(0, k) = B(2, k + 1) = dN(j, 0);
326 dNv(k + 1) = B(1, k + 1) = B(2, k) = dN(j, 1);
330 dNv(6) = B(0, 6) = B(2, 7) = 27. * ( dN(0, 0) *
N(1) *
N(2) +
N(0) * dN(1, 0) *
N(2) +
N(0) *
N(1) * dN(2, 0) );
331 dNv(7) = B(1, 7) = B(2, 6) = 27. * ( dN(0, 1) *
N(1) *
N(2) +
N(0) * dN(1, 1) *
N(2) +
N(0) *
N(1) * dN(2, 1) );
407 for (
int i = 1; i <= n.
giveSize(); i++ ) {
408 answer(0) += n.
at(i) * velocities.
at(i*2-1);
409 answer(1) += n.
at(i) * velocities.
at(i*2);
411 answer(0) += n.
at(1) * n.
at(2) * n.
at(3) * velocities.
at(7);
412 answer(1) += n.
at(1) * n.
at(2) * n.
at(3) * velocities.
at(8);
414 for (
int i = 1; i <= n_lin.
giveSize(); i++ ) {
415 answer(2) += n_lin.
at(i) * pressures.
at(i);
CrossSection * giveCrossSection()
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.
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 evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
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.
static FEI2dTrLin interp
Interpolation for pressure.
Abstract base class for all fluid materials.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
std::unique_ptr< ElementDofManager > bubble
The extra dofs from the bubble.
Class implementing internal element dof manager having some DOFs.
The element interface required by ZZNodalRecoveryModel.
bcGeomType
Type representing the geometric character of loading.
double & at(int i)
Coefficient access function.
virtual ~Tr1BubbleStokes()
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
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)
virtual FEInterpolation * giveInterpolation() const
void clear()
Clears receiver (zero size).
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Base class for dof managers.
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 updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Class representing a general abstraction for finite element interpolation class.
Class representing "master" degree of freedom.
virtual void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
static IntArray conservation_ordering
static IntArray edge_ordering[3]
Ordering of dofs on edges. Used to assemble edge loads.
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 Interface * giveInterface(InterfaceType it)
Interface requesting service.
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 giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep)
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).
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation Jacobian on the requested boundary.
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 int giveApproxOrder()=0
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
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 .
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
This abstract class represent a general base element class for fluid dynamic problems.
virtual void computeField(ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
Implementation of matrix containing floating point numbers.
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
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.
static IntArray momentum_ordering
Ordering of dofs in element. Used to assemble the element stiffness.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
void add(const FloatMatrix &a)
Adds matrix to the receiver.
void zero()
Zeroes all coefficients of receiver.
virtual void giveInternalDofManDofIDMask(int i, IntArray &answer) const
Returns internal dofmanager dof mask for node.
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual int giveNumberOfInternalDofManagers() const
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...
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
void zero()
Zeroes all coefficient of receiver.
InterfaceType
Enumerative type, used to identify interface type.
Tr1BubbleStokes(int n, Domain *d)
virtual DofManager * giveInternalDofManager(int i) const
Returns i-th internal element dof manager of the receiver.
Load is base abstract class for all loads.
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
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 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)...
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 double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Class representing integration point in finite element program.
IntArray boundaryLoadArray
Class representing solution step.
int numberOfDofMans
Number of dofmanagers.
void add(const FloatArray &src)
Adds array src to receiver.
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .