OOFEM 3.0
Loading...
Searching...
No Matches
surfacetensionbc.C
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 - 2025 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#include "surfacetensionbc.h"
36#include "element.h"
37#include "node.h"
38#include "crosssection.h"
39#include "error.h"
40#include "feinterpol.h"
41#include "feinterpol2d.h"
42#include "feinterpol3d.h"
43#include "classfactory.h"
44#include "set.h"
45#include "sparsemtrx.h"
46
47#include "integrationdomain.h"
48#include "integrationrule.h"
49#include "gausspoint.h"
50#include "mathfem.h"
51
52#include <utility>
53#include <list>
54#include <memory>
55
56#ifdef _OPENMP
57#include <omp.h>
58#endif
59
60namespace oofem {
62
63void SurfaceTensionBoundaryCondition :: initializeFrom(InputRecord &ir)
64{
65 ActiveBoundaryCondition :: initializeFrom(ir);
66
68
70}
71
72void SurfaceTensionBoundaryCondition :: giveLocationArrays(std :: vector< IntArray > &rows, std :: vector< IntArray > &cols, CharType type,
74{
75 if ( !this->useTangent || type != TangentStiffnessMatrix ) {
76 return;
77 }
78
79 Set *set = this->giveDomain()->giveSet(this->set);
80 const IntArray &boundaries = set->giveBoundaryList();
81
82 rows.resize(boundaries.giveSize() / 2);
83 cols.resize(boundaries.giveSize() / 2);
84
85 for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
86 Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
87 int boundary = boundaries.at(pos * 2);
88
89 const auto &bNodes = e->giveInterpolation()->boundaryGiveNodes(boundary, e->giveGeometryType());
90
91 e->giveBoundaryLocationArray(rows [ pos ], bNodes, this->dofs, r_s);
92 e->giveBoundaryLocationArray(cols [ pos ], bNodes, this->dofs, c_s);
93 }
94}
95
96void SurfaceTensionBoundaryCondition :: assemble(SparseMtrx &answer, TimeStep *tStep,
97 CharType type, const UnknownNumberingScheme &r_s,
98 const UnknownNumberingScheme &c_s, double scale,
99 void*lock)
100{
101 if ( !this->useTangent || type != TangentStiffnessMatrix ) {
102 return;
103 }
104
105 OOFEM_ERROR("Not implemented yet.");
106#if 0
107 FloatMatrix Ke;
108 IntArray r_loc, c_loc;
109
110 Set *set = this->giveDomain()->giveSet(this->set);
111 const IntArray &boundaries = set->giveBoundaryList();
112
113 for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
114 Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
115 int boundary = boundaries.at(pos * 2);
116
117 const auto &bNodes = e->giveInterpolation()->boundaryGiveNodes(boundary, e->giveGeometryType());
118
119 e->giveBoundaryLocationArray(r_loc, bNodes, this->dofs, r_s);
120 e->giveBoundaryLocationArray(c_loc, bNodes, this->dofs, c_s);
121 this->computeTangentFromElement(Ke, e, boundary, tStep);
122 Ke.times(scale);
123#ifdef _OPENMP
124 if (lock) omp_set_lock(static_cast<omp_lock_t*>(lock));
125#endif
126 answer.assemble(r_loc, c_loc, Ke);
127#ifdef _OPENMP
128 if (lock) omp_unset_lock(static_cast<omp_lock_t*>(lock));
129#endif
130 }
131#endif
132}
133
134void SurfaceTensionBoundaryCondition :: assembleVector(FloatArray &answer, TimeStep *tStep,
135 CharType type, ValueModeType mode,
136 const UnknownNumberingScheme &s,
137 FloatArray *eNorms,
138 void*lock)
139{
140 if ( type != ExternalForcesVector ) {
141 return;
142 }
143
144 FloatArray fe;
145 IntArray loc, masterdofids;
146
147 Set *set = this->giveDomain()->giveSet(this->set);
148 const IntArray &boundaries = set->giveBoundaryList();
149
150 for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
151 Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
152 int boundary = boundaries.at(pos * 2);
153
154 const auto &bNodes = e->giveInterpolation()->boundaryGiveNodes(boundary, e->giveGeometryType());
155
156 e->giveBoundaryLocationArray(loc, bNodes, this->dofs, s, & masterdofids);
157 this->computeLoadVectorFromElement(fe, e, boundary, tStep);
158#ifdef _OPENMP
159 if (lock) omp_set_lock(static_cast<omp_lock_t*>(lock));
160#endif
161 answer.assemble(fe, loc);
162 if ( eNorms ) {
163 eNorms->assembleSquared(fe, masterdofids);
164 }
165#ifdef _OPENMP
166 if (lock) omp_unset_lock(static_cast<omp_lock_t*>(lock));
167#endif
168 }
169}
170
171void SurfaceTensionBoundaryCondition :: computeTangentFromElement(FloatMatrix &answer, Element *e, int side, TimeStep *tStep)
172{
174 if ( !fei ) {
175 OOFEM_ERROR("No interpolation available for element.");
176 }
177 std :: unique_ptr< IntegrationRule > iRule( fei->giveBoundaryIntegrationRule(fei->giveInterpolationOrder()-1, side, e->giveGeometryType()) );
178
180 int nodes = e->giveNumberOfNodes();
181 if ( side == -1 ) {
182 side = 1;
183 }
184
185 answer.clear();
186
187 if ( nsd == 2 ) {
188 FEInterpolation2d *fei2d = static_cast< FEInterpolation2d * >(fei);
189
191 FloatMatrix xy(2, nodes);
192 Node *node;
193 for ( int i = 1; i <= nodes; i++ ) {
194 node = e->giveNode(i);
195 xy.at(1, i) = node->giveCoordinate(1);
196 xy.at(2, i) = node->giveCoordinate(2);
197 }
198
199 FloatArray tmpA(2 *nodes);
200 FloatArray es; // Tangent vector to curve
201 FloatArray dNds;
202 FloatMatrix B(2, 2 *nodes);
203 B.zero();
204
205 if ( e->giveDomain()->isAxisymmetric() ) {
207 FloatArray gcoords;
208 FloatArray tmpB(2 *nodes);
209 for ( GaussPoint *gp: *iRule ) {
210 fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
211 fei->boundaryEvalN( N, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
212 double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
213 fei->boundaryLocal2Global( gcoords, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
214 double r = gcoords(0); // First coordinate is the radial coord.
215
216 es.beProductOf(xy, dNds);
217
218 // Construct the different matrices in the integrand;
219 for ( int i = 0; i < nodes; i++ ) {
220 tmpA(i * 2 + 0) = dNds(i) * es(0);
221 tmpA(i * 2 + 1) = dNds(i) * es(1);
222 tmpB(i * 2 + 0) = N(i);
223 tmpB(i * 2 + 1) = 0;
224 B(i * 2, 0) = B(i * 2 + 1, 1) = dNds(i);
225 }
226
227 double dV = 2 *M_PI *gamma *J *gp->giveWeight();
228 answer.plusDyadUnsym(tmpA, tmpB, dV);
229 answer.plusDyadUnsym(tmpB, tmpA, dV);
230 answer.plusProductSymmUpper(B, B, r * dV);
231 answer.plusDyadUnsym(tmpA, tmpA, -r * dV);
232 }
233 } else {
234 for ( GaussPoint *gp: *iRule ) {
235 double t = e->giveCrossSection()->give(CS_Thickness, gp);
236 fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
237 double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
238
239 es.beProductOf(xy, dNds);
240
241 // Construct the different matrices in the integrand;
242 for ( int i = 0; i < nodes; i++ ) {
243 tmpA(i * 2 + 0) = dNds(i) * es(0);
244 tmpA(i * 2 + 1) = dNds(i) * es(1);
245 B(i * 2, 0) = B(i * 2 + 1, 1) = dNds(i);
246 }
247
248 double dV = t * gamma * J * gp->giveWeight();
249 answer.plusProductSymmUpper(B, B, dV);
250 answer.plusDyadSymmUpper(tmpA, -dV);
251 }
252 }
253
254 answer.symmetrized();
255 } else if ( nsd == 3 ) {
256
257 //FEInterpolation3d *fei3d = static_cast< FEInterpolation3d * >(fei);
258
259 OOFEM_ERROR("3D tangents not implemented yet.");
260#if 0
261 //FloatMatrix tmp(3 *nodes, 3 *nodes);
262 FloatMatrix dNdx;
263 FloatArray n;
264 for ( GaussPoint *gp: *iRule ) {
265 fei3d->surfaceEvaldNdx( dNdx, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
266 /*double J = */ fei->boundaryEvalNormal( n, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
267 //double dV = gamma * J * gp->giveWeight();
268
269 for ( int i = 0; i < nodes; i++ ) {
270 //tmp(3*i+0) = dNdx(i,0) - (dNdx(i,0)*n(0)* + dNdx(i,1)*n(1) + dNdx(i,2)*n(2))*n(0);
271 //tmp(3*i+1) = dNdx(i,1) - (dNdx(i,0)*n(0)* + dNdx(i,1)*n(1) + dNdx(i,2)*n(2))*n(1);
272 //tmp(3*i+2) = dNdx(i,2) - (dNdx(i,0)*n(0)* + dNdx(i,1)*n(1) + dNdx(i,2)*n(2))*n(2);
273 }
274 //answer.plusProductSymmUpper(A,B, dV);
276 }
277#endif
278 } else {
279 OOFEM_WARNING("Only 2D or 3D is possible!");
280 }
281}
282
283void SurfaceTensionBoundaryCondition :: computeLoadVectorFromElement(FloatArray &answer, Element *e, int side, TimeStep *tStep)
284{
286 if ( !fei ) {
287 OOFEM_ERROR("No interpolation or default integration available for element.");
288 }
289 std :: unique_ptr< IntegrationRule > iRule( fei->giveBoundaryIntegrationRule(fei->giveInterpolationOrder()-1, side, e->giveGeometryType()) );
290
292
293 if ( side == -1 ) {
294 side = 1;
295 }
296
297 answer.clear();
298
299 if ( nsd == 2 ) {
300
301 FEInterpolation2d *fei2d = static_cast< FEInterpolation2d * >(fei);
302
304 const auto &bnodes = fei2d->boundaryGiveNodes(side, e->giveGeometryType());
305 int nodes = bnodes.giveSize();
306 FloatMatrix xy(2, nodes);
307 for ( int i = 1; i <= nodes; i++ ) {
308 Node *node = e->giveNode(bnodes.at(i));
310 xy.at(1, i) = node->giveCoordinate(1);
311 xy.at(2, i) = node->giveCoordinate(2);
312 }
313
314 FloatArray tmp; // Integrand
315 FloatArray es; // Tangent vector to curve
316 FloatArray dNds;
317
318 if ( e->giveDomain()->isAxisymmetric() ) {
320 FloatArray gcoords;
321 for ( GaussPoint *gp: *iRule ) {
322 fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
323 fei->boundaryEvalN( N, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
324 double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
325 fei->boundaryLocal2Global( gcoords, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
326 double r = gcoords(0); // First coordinate is the radial coord.
327
328 es.beProductOf(xy, dNds);
329
330 tmp.resize( 2 * nodes);
331 for ( int i = 0; i < nodes; i++ ) {
332 tmp(2 * i) = dNds(i) * es(0) * r + N(i);
333 tmp(2 * i + 1) = dNds(i) * es(1) * r;
334 }
335
336 answer.add(- 2 * M_PI * gamma * J * gp->giveWeight(), tmp);
337 }
338 } else {
339 for ( GaussPoint *gp: *iRule ) {
340 double t = e->giveCrossSection()->give(CS_Thickness, gp);
341 fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
342 double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
343 es.beProductOf(xy, dNds);
344
345 tmp.resize( 2 * nodes);
346 for ( int i = 0; i < nodes; i++ ) {
347 tmp(2 * i) = dNds(i) * es(0);
348 tmp(2 * i + 1) = dNds(i) * es(1);
349 //B.at(1, 1+i*2) = B.at(2, 2+i*2) = dNds(i);
350 }
351 //tmp.beTProductOf(B, es);
352
353 answer.add(- t * gamma * J * gp->giveWeight(), tmp);
354 }
355 }
356 } else if ( nsd == 3 ) {
357
358 FEInterpolation3d *fei3d = static_cast< FEInterpolation3d * >(fei);
359 FloatArray n, surfProj;
360 FloatMatrix dNdx, B;
361 for ( GaussPoint *gp: *iRule ) {
362 fei3d->surfaceEvaldNdx( dNdx, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
363 double J = fei->boundaryEvalNormal( n, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
364
365 // [I - n(x)n] in voigt form:
366 surfProj = Vec6(1. - n(0)*n(0), 1. - n(1)*n(1), 1. - n(2)*n(2),
367 - n(1)*n(2), - n(0)*n(2), - n(0)*n(1)
368 );
369
370 // Construct B matrix of the surface nodes
371 B.resize(6, 3 * dNdx.giveNumberOfRows());
372 for ( int i = 1; i <= dNdx.giveNumberOfRows(); i++ ) {
373 B.at(1, 3 * i - 2) = dNdx.at(i, 1);
374 B.at(2, 3 * i - 1) = dNdx.at(i, 2);
375 B.at(3, 3 * i - 0) = dNdx.at(i, 3);
376
377 B.at(5, 3 * i - 2) = B.at(4, 3 * i - 1) = dNdx.at(i, 3);
378 B.at(6, 3 * i - 2) = B.at(4, 3 * i - 0) = dNdx.at(i, 2);
379 B.at(6, 3 * i - 1) = B.at(5, 3 * i - 0) = dNdx.at(i, 1);
380 }
381
382 answer.plusProduct(B, surfProj, -gamma * J * gp->giveWeight() );
383 }
384 }
385}
386} // end namespace oofem
#define N(a, b)
#define REGISTER_BoundaryCondition(class)
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
double giveCoordinate(int i) const
Definition dofmanager.h:383
bool isAxisymmetric()
Returns true of axisymmetry is in effect.
Definition domain.C:1144
int giveNumberOfSpatialDimensions()
Returns number of spatial dimensions.
Definition domain.C:1137
Node * giveNode(int i) const
Definition element.h:629
virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
Definition element.C:485
virtual FEInterpolation * giveInterpolation() const
Definition element.h:648
virtual int giveNumberOfNodes() const
Definition element.h:703
CrossSection * giveCrossSection()
Definition element.C:534
virtual Element_Geometry_Type giveGeometryType() const =0
virtual void edgeEvaldNds(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
IntArray boundaryGiveNodes(int boundary, const Element_Geometry_Type) const override
virtual void surfaceEvaldNdx(FloatMatrix &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
virtual std::unique_ptr< IntegrationRule > giveBoundaryIntegrationRule(int order, int boundary, const Element_Geometry_Type) const
Definition feinterpol.C:101
virtual IntArray boundaryGiveNodes(int boundary, const Element_Geometry_Type) const =0
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
int giveInterpolationOrder() const
Definition feinterpol.h:199
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double boundaryEvalNormal(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void boundaryEvalN(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
Domain * giveDomain() const
Definition femcmpnn.h:97
void assemble(const FloatArray &fe, const IntArray &loc)
Definition floatarray.C:616
void resize(Index s)
Definition floatarray.C:94
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Definition floatarray.C:288
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Definition floatarray.C:637
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double f)
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
void zero()
Zeroes all coefficient of receiver.
void plusDyadSymmUpper(const FloatArray &a, double dV)
int giveNumberOfRows() const
Returns number of rows of receiver.
double at(std::size_t i, std::size_t j) const
int set
Set number for boundary condition to be applied to.
IntArray dofs
Dofs that b.c. is applied to (relevant for Dirichlet type b.c.s).
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
bool useTangent
Determines if tangent should be used.
void computeLoadVectorFromElement(FloatArray &answer, Element *e, int side, TimeStep *tStep)
void computeTangentFromElement(FloatMatrix &answer, Element *e, int side, TimeStep *tStep)
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
#define M_PI
Definition mathfem.h:52
@ CS_Thickness
Thickness.
static FloatArray Vec6(const double &a, const double &b, const double &c, const double &d, const double &e, const double &f)
Definition floatarray.h:610
#define _IFT_SurfaceTensionBoundaryCondition_useTangent
#define _IFT_SurfaceTensionBoundaryCondition_gamma

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011