35 #ifndef tr1_2d_supg2_h 36 #define tr1_2d_supg2_h 45 #define _IFT_TR1_2D_SUPG2_Name "tr1supg2" 141 const FloatArray &normal,
const double p,
bool updFlag);
179 const FloatArray &normal,
const double p,
bool updFlag);
185 #endif // tr1_2d_supg2_h virtual void giveElementCenter(LEPlic *mat_interface, FloatArray ¢er, bool updFlag)
Computes the receiver center (in updated Lagrangian configuration).
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual Element * giveElement()
Return number of receiver's element.
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual const char * giveInputRecordName() const
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
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...
virtual void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Abstract class representing field of primary variables (those, which are unknown and are typically as...
Material * _giveMaterial(int indx)
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes SLIC stabilization term for momentum balance equation(s).
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
Computes corresponding volume fraction to given interface position.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
Class representing 2d linear triangular element for solving incompressible fluid with SUPG solver...
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Abstract base class for all finite elements.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Class implementing an array of integers.
MatResponseMode
Describes the character of characteristic material matrix.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Abstract base class representing integration rule.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
int mat[2]
mat[0] reference fluid.
void updateIntegrationRules()
virtual int giveDefaultIntegrationRule() const
Returns id of default integration rule.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Material * giveMaterial(int n)
Service for accessing particular domain material model.
virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for momentum balance equations(s).
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the linear advection term for mass conservation equation.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
TR1_2D_SUPG2(int n, Domain *d)
virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag)
Assembles receiver volume.
virtual double computeMyVolume(LEPlic *matInterface, bool updFlag)
Computes the volume of receiver.
virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for momentum balance equations(s) with respect to nodal ve...
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...
#define _IFT_TR1_2D_SUPG2_Name
Abstract base class for all material models.
virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for mass conservation equation with respect to nodal veloc...
virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms (generalized mass matrix with stabilization terms) for momentum balance e...
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.
virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms for mass conservation equation.
double computeVolumeAroundID(GaussPoint *gp, integrationDomain id, const std::vector< FloatArray > &idpoly)
Abstract base class representing Lagrangian-Eulerian (moving) material interfaces.
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).
virtual const char * giveClassName() const
Class representing vector of real numbers.
Class representing 2D polygon.
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *)
Implementation of matrix containing floating point numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for mass conservation equation.
virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume)
Truncates given material polygon to receiver.
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
virtual void updateElementForNewInterfacePosition(TimeStep *tStep)
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
void zero()
Zeroes all coefficients of receiver.
Polygon myPoly[2]
myPoly[0] occupied by reference fluid.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for momentum balance equations(s).
std::vector< FloatArray > vcoords[2]
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) ...
void zero()
Zeroes all coefficient of receiver.
InterfaceType
Enumerative type, used to identify interface type.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Class representing 2d linear triangular element for solving incompressible fluid with SUPG solver...
virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes diffusion derivative terms for mass conservation equation.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void postInitialize()
Performs post initialization steps.
double p
Line constant of line segment representing interface.
std::unique_ptr< IntegrationRule > defaultIRule
default integration rule over element volume.
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
void computeNVector(FloatArray &answer, GaussPoint *gp)
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Updates the stabilization coefficients used for CBS and SUPG algorithms.
Class representing integration point in finite element program.
FloatArray normal
Interface segment normal.
Class representing solution step.
InternalStateMode
Determines the mode of internal variable.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
virtual double computeCriticalTimeStep(TimeStep *tStep)
Computes the critical time increment.
virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes advection terms for mass conservation equation.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for mass conservation equation.
virtual void updateVolumePolygons(Polygon &referenceFluidPoly, Polygon &secondFluidPoly, int &rfPoints, int &sfPoints, const FloatArray &normal, const double p, bool updFlag)
virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes nonlinear advection terms for momentum balance equations(s).
void resize(int s)
Resizes receiver towards requested size.