73 if ( !this->
useTangent || type != TangentStiffnessMatrix ) {
78 const IntArray &boundaries =
set->giveBoundaryList();
82 cols.resize(boundaries.
giveSize() / 2);
84 for (
int pos = 1; pos <= boundaries.
giveSize() / 2; ++pos ) {
86 int boundary = boundaries.
at(pos * 2);
98 if ( !this->
useTangent || type != TangentStiffnessMatrix ) {
108 const IntArray &boundaries =
set->giveBoundaryList();
110 for (
int pos = 1; pos <= boundaries.
giveSize() / 2; ++pos ) {
112 int boundary = boundaries.
at(pos * 2);
128 if ( type != ExternalForcesVector ) {
136 const IntArray &boundaries =
set->giveBoundaryList();
138 for (
int pos = 1; pos <= boundaries.
giveSize() / 2; ++pos ) {
140 int boundary = boundaries.
at(pos * 2);
157 OOFEM_ERROR(
"No interpolation available for element.");
175 for (
int i = 1; i <= nodes; i++ ) {
196 double r = gcoords(0);
201 for (
int i = 0; i < nodes; i++ ) {
202 tmpA(i * 2 + 0) = dNds(i) * es(0);
203 tmpA(i * 2 + 1) = dNds(i) * es(1);
204 tmpB(i * 2 + 0) =
N(i);
206 B(i * 2, 0) = B(i * 2 + 1, 1) = dNds(i);
209 double dV = 2 *
M_PI *
gamma *J *gp->giveWeight();
224 for (
int i = 0; i < nodes; i++ ) {
225 tmpA(i * 2 + 0) = dNds(i) * es(0);
226 tmpA(i * 2 + 1) = dNds(i) * es(1);
227 B(i * 2, 0) = B(i * 2 + 1, 1) = dNds(i);
230 double dV = t *
gamma * J * gp->giveWeight();
237 }
else if ( nsd == 3 ) {
251 for (
int i = 0; i < nodes; i++ ) {
268 OOFEM_ERROR(
"No interpolation or default integration available for element.");
289 for (
int i = 1; i <= nodes; i++ ) {
308 double r = gcoords(0);
313 for (
int i = 0; i < nodes; i++ ) {
314 tmp(2 * i) = dNds(i) * es(0) * r +
N(i);
315 tmp(2 * i + 1) = dNds(i) * es(1) * r;
318 answer.
add(- 2 *
M_PI *
gamma * J * gp->giveWeight(), tmp);
328 for (
int i = 0; i < nodes; i++ ) {
329 tmp(2 * i) = dNds(i) * es(0);
330 tmp(2 * i + 1) = dNds(i) * es(1);
335 answer.
add(- t *
gamma * J * gp->giveWeight(), tmp);
338 }
else if ( nsd == 3 ) {
348 surfProj = {1. - n(0)*n(0), 1. - n(1)*n(1), 1. - n(2)*n(2),
349 - n(1)*n(2), - n(0)*n(2), - n(0)*n(1),
355 B.
at(1, 3 * i - 2) = dNdx.
at(i, 1);
356 B.
at(2, 3 * i - 1) = dNdx.
at(i, 2);
357 B.
at(3, 3 * i - 0) = dNdx.
at(i, 3);
359 B.
at(5, 3 * i - 2) = B.
at(4, 3 * i - 1) = dNdx.
at(i, 3);
360 B.
at(6, 3 * i - 2) = B.
at(4, 3 * i - 0) = dNdx.
at(i, 2);
361 B.
at(6, 3 * i - 1) = B.
at(5, 3 * i - 0) = dNdx.
at(i, 1);
CrossSection * giveCrossSection()
#define _IFT_SurfaceTensionBoundaryCondition_gamma
virtual void assembleVector(FloatArray &answer, TimeStep *tStep, CharType type, ValueModeType mode, const UnknownNumberingScheme &s, FloatArray *eNorms=NULL)
Assembles B.C.
REGISTER_BoundaryCondition(BoundaryCondition)
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
Assembles sparse matrix from contribution of local elements.
Base class for all matrices stored in sparse format.
IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void boundaryGiveNodes(IntArray &answer, int boundary)
Gives the boundary nodes for requested boundary number.
virtual double boundaryEvalNormal(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the normal on the requested boundary.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Class representing a general abstraction for surface finite element interpolation class...
void clear()
Clears receiver (zero size).
int giveInterpolationOrder()
Returns the interpolation order.
bool isAxisymmetric()
Returns true of axisymmetry is in effect.
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale=1.0)
Assembles B.C.
Abstract base class for all finite elements.
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
Adds to the receiver the product .
virtual double giveCoordinate(int i)
int giveNumberOfSpatialDimensions()
Returns number of spatial dimensions.
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual FEInterpolation * giveInterpolation() const
void computeLoadVectorFromElement(FloatArray &answer, Element *e, int side, TimeStep *tStep)
Helper function for computing the contributions to the load vector.
virtual int giveNumberOfNodes() const
Returns number of nodes of receiver.
Class representing a general abstraction for finite element interpolation class.
Element * giveElement(int n)
Service for accessing particular domain fe element.
double gamma
Surface tension.
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Maps the local boundary coordinates to global.
void times(double f)
Multiplies receiver by factor f.
Set of elements, boundaries, edges and/or nodes.
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
Wrapper around element definition to provide FEICellGeometry interface.
Set * giveSet(int n)
Service for accessing particular domain set.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the determinant of the transformation Jacobian on the requested boundary.
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
IntArray dofs
Dofs that b.c. is applied to (relevant for Dirichlet type b.c.s).
double at(int i, int j) const
Coefficient access function.
void resize(int n)
Checks size of receiver towards requested bounds.
virtual void scale(double s)
Scales the receiver according to given value.
virtual void boundaryGiveNodes(IntArray &answer, int boundary)=0
Gives the boundary nodes for requested boundary number.
Class representing vector of real numbers.
void plusDyadSymmUpper(const FloatArray &a, double dV)
Adds to the receiver the dyadic product .
Implementation of matrix containing floating point numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
virtual void surfaceEvaldNdx(FloatMatrix &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point...
virtual void giveLocationArrays(std::vector< IntArray > &rows, std::vector< IntArray > &cols, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
Gives a list of location arrays that will be assembled.
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.
Class representing a general abstraction for surface finite element interpolation class...
virtual IntegrationRule * giveBoundaryIntegrationRule(int order, int boundary)
Sets up a suitable integration rule for integrating over the requested boundary.
void computeTangentFromElement(FloatMatrix &answer, Element *e, int side, TimeStep *tStep)
Helper function for computing the tangent ( )
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Assembles the array fe with each component squared.
virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
Returns the location array for the boundary of the element.
virtual void edgeEvaldNds(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point...
void zero()
Zeroes all coefficient of receiver.
Domain * giveDomain() const
#define _IFT_SurfaceTensionBoundaryCondition_useTangent
the oofem namespace is to define a context or scope in which all oofem names are defined.
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Class implementing node in finite element mesh.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
bool useTangent
Determines if tangent should be used.
int giveNumberOfRows() const
Returns number of rows of receiver.
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
virtual void boundaryEvalN(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the basis functions on the requested boundary.
Class representing integration point in finite element program.
#define OOFEM_WARNING(...)
Class representing solution step.
void add(const FloatArray &src)
Adds array src to receiver.
void resize(int s)
Resizes receiver towards requested size.
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .