OOFEM 3.0
Loading...
Searching...
No Matches
structuralinterfaceelementphf.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
38
40#include "feinterpol.h"
41#include "domain.h"
42#include "material.h"
43#include "gausspoint.h"
45#include "intarray.h"
46#include "floatarray.h"
47#include "floatmatrix.h"
48#include "mathfem.h"
49#include "dof.h"
50
51
52namespace oofem {
53StructuralInterfaceElementPhF :: StructuralInterfaceElementPhF(int n, Domain *aDomain) : StructuralInterfaceElement(n, aDomain)
54{
55 this->internalLength = 1.0e-1;
56}
57
58
59// remove later
60void
61StructuralInterfaceElementPhF :: computeDisplacementUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *tStep)
62{
63 IntArray dofIdArray;
64 this->giveDofManDofIDMask_u(dofIdArray);
65 this->computeVectorOf(dofIdArray, valueMode, tStep, answer);
66}
67
68
69void
70StructuralInterfaceElementPhF :: computeDamageUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *tStep)
71{
72 IntArray dofIdArray;
73 this->giveDofManDofIDMask_d(dofIdArray);
74 this->computeVectorOf(dofIdArray, valueMode, tStep, answer);
75}
76
77
78void
79StructuralInterfaceElementPhF :: computeLocationArrayOfDofIDs( const IntArray &dofIdArray, IntArray &answer )
80{
81 // Routine to compute the local ordering array an element given a dofid array.
82 answer.resize( 0 );
83 int k = 0;
84 for(int i = 1; i <= this->giveNumberOfDofManagers(); i++) {
85 DofManager *dMan = this->giveDofManager( i );
86 for( int j = 1; j <= dofIdArray.giveSize( ); j++ ) {
87
88 if ( dMan->hasDofID( (DofIDItem) dofIdArray.at( j ) ) ) {
89 // hack
90 answer.followedBy( k + j );
91 }
92 }
93 k += dMan->giveNumberOfDofs( );
94 }
95}
96
97
98void
99StructuralInterfaceElementPhF :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
100{
101 // Computes the internal forces corresponding to the two fields u & d
102 IntArray IdMask_u, IdMask_d;
103 this->giveDofManDofIDMask_u( IdMask_u );
104 this->giveDofManDofIDMask_d( IdMask_d );
107 int ndofs = this->computeNumberOfDofs();
108 answer.resize( ndofs);
109 answer.zero();
110 FloatArray answer_u(0);
111 FloatArray answer_d(0);
112 this->giveInternalForcesVector_u(answer_u, tStep, useUpdatedGpRecord);
113 this->giveInternalForcesVector_d(answer_d, tStep, useUpdatedGpRecord);
114 answer.assemble(answer_u, loc_u);
115 answer.assemble(answer_d, loc_d);
116}
117
118
119void
120StructuralInterfaceElementPhF :: giveInternalForcesVector_u(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
121{
123 FloatArray u, traction, jump;
124
125 IntArray dofIdArray;
126 this->giveDofManDofIDMask_u(dofIdArray);
127 this->computeVectorOf(dofIdArray, VM_Total, tStep, u);
128
129 // subtract initial displacements, if defined
130 if ( initialDisplacements.giveSize() ) {
132 }
133
134 // zero answer will resize accordingly when adding first contribution
135 answer.clear();
136
137 for ( auto &ip: *this->giveDefaultIntegrationRulePtr() ) {
138 this->computeNmatrixAt(ip, N);
139
140 jump.beProductOf(N, u);
141 this->computeTraction(traction, ip, jump, tStep);
142 //traction.resize(2);
143
144 // compute internal cohesive forces as f = N^T*traction dA
145 double dA = this->computeAreaAround(ip);
146 answer.plusProduct(N, traction, dA);
147 }
148}
149
150void
151StructuralInterfaceElementPhF :: giveInternalForcesVector_d(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
152{
153 // computes int_A ( N^t * ( d + l*f ) + B^t * l^2 * gradd(d) )* dA
154 FloatArray BStress, Nd, BS, a_d, grad_d;
155 FloatMatrix B;
156 answer.clear();
157 this->computeDamageUnknowns( a_d, VM_Total, tStep );
158
159 double l = this->giveInternalLength();
160
161 for ( auto &ip: *this->giveDefaultIntegrationRulePtr() ) {
162
163 // Part 1
164 StructuralInterfaceMaterialPhF *mat = static_cast< StructuralInterfaceMaterialPhF * >( this->giveInterfaceCrossSection()->giveInterfaceMaterial() );
165
166 double d = computeDamageAt( ip, VM_Total, tStep );
167 double facN = d + l * mat->giveDrivingForce(ip);
168
169 this->giveInterpolation()->evalN(Nd, ip->giveNaturalCoordinates(), FEIElementGeometryWrapper(this));
170
171 double dA = this->computeAreaAround(ip);
172 answer.add(facN*dA, Nd);
173
174 // Part 2
175 this->computeBd_matrixAt(ip, B );
176 grad_d.beProductOf(B, { a_d.at(1), a_d.at(2) });
177 BStress = grad_d * l * l;
178 BS.beTProductOf(B, BStress);
179 answer.add(dA, BS);
180 }
181
182 auto temp = answer;
183 answer.resize(4);
184 answer.zero();
185 answer.assemble(temp,{1,2});
186 answer.assemble(temp,{3,4});
187}
188
189
190double
191StructuralInterfaceElementPhF :: computeDamageAt(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
192{
193 // d = N_d * a_d
194 FloatArray dVec;
195 computeDamageUnknowns(dVec, valueMode, stepN);
196 FloatArray Nvec;
197 this->giveInterpolation()->evalN(Nvec, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this));
198 //return Nvec.dotProduct(dVec);
199 return Nvec.dotProduct( {dVec.at(1), dVec.at(2) });
200}
201
202
203void
204StructuralInterfaceElementPhF :: computeNd_matrixAt(const FloatArray &lCoords, FloatMatrix &N)
205{
206 FloatArray Nvec;
207 this->giveInterpolation( )->evalN( Nvec, lCoords, FEIElementGeometryWrapper( this ) );
208 //N.resize(1, Nvec.giveSize());
209 N.beNMatrixOf(Nvec,1);
210}
211
212void
213StructuralInterfaceElementPhF :: computeBd_matrixAt(GaussPoint *gp, FloatMatrix &answer)
214{
215 FloatMatrix dNdxi;
216 FloatMatrix G;
217 this->computeCovarBaseVectorsAt(gp, G);
218
219 this->giveInterpolation( )->evaldNdxi( dNdxi, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper( this ) );
220 answer.beProductTOf(G,dNdxi);
221 //answer.beTranspositionOf( dNdx );
222}
223
224
225int
226StructuralInterfaceElementPhF :: computeNumberOfDofs()
227{
228 //NLStructuralElement *el = static_cast< NLStructuralElement* >(this->giveElement( ) );
229 int nDofs = 0;
230 for( int i = 1; i <= this->giveNumberOfDofManagers(); i++) {
231 nDofs += this->giveDofManager( i )->giveNumberOfDofs();
232 }
233 return nDofs;
234}
235
236
237void
238StructuralInterfaceElementPhF :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
239{
240 //set displacement and nonlocal location array
242 //IntArray IdMask_u, IdMask_d;
243 //this->giveDofManDofIDMask_u( IdMask_u );
244 //this->giveDofManDofIDMask_d( IdMask_d );
245
248
249 int nDofs = this->computeNumberOfDofs();
250 answer.resize( nDofs, nDofs );
251 answer.zero();
252
253 FloatMatrix answer1, answer2, answer3, answer4;
254
255 this->computeStiffnessMatrix_uu(answer1, rMode, tStep);
256 //StructuralInterfaceElement :: computeStiffnessMatrix(answer1, rMode, tStep);
257
258 //this->computeStiffnessMatrix_ud(answer2, rMode, tStep);
259 //this->computeStiffnessMatrix_du(answer3, rMode, tStep); //symmetric
260 //answer3.beTranspositionOf( answer2 );
261
262 this->computeStiffnessMatrix_dd(answer4, rMode, tStep);
263
264 answer.assemble( answer1, loc_u, loc_u );
265 //answer.assemble( answer2, loc_u, loc_d );
266 //answer.assemble( answer3, loc_d, loc_u );
267 answer.assemble( answer4, loc_d, loc_d );
268}
269
270
271void
272StructuralInterfaceElementPhF :: computeStiffnessMatrix_uu(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
273{
274 // This is the regular stiffness matrix times G
275
276 // Computes the stiffness matrix of the receiver K_cohesive = int_A ( N^t * D * N ) dA
277 // with the stiffness D = dT/dj
278 FloatMatrix N, D, DN;
279 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
280 answer.clear();
281
282 FloatMatrix rotationMatGtoL;
283 for ( auto &ip: *this->giveDefaultIntegrationRulePtr() ) {
284
285 this->giveStiffnessMatrix_Eng(D, rMode, ip, tStep);
286
287 this->computeTransformationMatrixAt(ip, rotationMatGtoL);
288 D.rotatedWith(rotationMatGtoL, 't'); // transform stiffness to global coord system
289
290 this->computeNmatrixAt(ip, N);
291 DN.beProductOf(D, N);
292 double dA = this->computeAreaAround(ip);
293 if ( matStiffSymmFlag ) {
294 answer.plusProductSymmUpper(N, DN, dA);
295 } else {
296 answer.plusProductUnsym(N, DN, dA);
297 }
298 }
299
300 if ( matStiffSymmFlag ) {
301 answer.symmetrized();
302 }
303}
304
305
306void
307StructuralInterfaceElementPhF :: computeStiffnessMatrix_ud(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
308{
309 FloatMatrix B, DB, N_d, DN, S(3,1);
310 FloatArray stress;
311#if 1
312 answer.clear();
313
314 for ( auto &ip: *this->giveDefaultIntegrationRulePtr() ) {
315 double dA = this->computeAreaAround(ip);
316
317 // compute int_V ( B^t * D_B * B )dV
318 this->computeNd_matrixAt(ip->giveNaturalCoordinates(), N_d);
319 // stress
320
321 FloatMatrix B_u, N;
322 this->computeNmatrixAt(ip, N);
323 //stress = matStat->giveTempStressVector();
324 //stress.times( this->computeGPrim( ip, VM_Total, tStep ) );
325 S.setColumn(stress,1);
326 DN.beProductOf(S, N_d);
327 answer.plusProductUnsym(N, DN, dA);
328 }
329#endif
330}
331
332
333void
334StructuralInterfaceElementPhF :: computeStiffnessMatrix_dd(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
335{
336 double l = this->giveInternalLength();
337
338 FloatMatrix B_d, N_d;
339 answer.clear();
340
341 FloatMatrix rotationMatGtoL;
342 for ( auto &ip: *this->giveDefaultIntegrationRulePtr() ) {
343 StructuralInterfaceMaterialPhF *mat = static_cast< StructuralInterfaceMaterialPhF * >( this->giveInterfaceCrossSection()->giveInterfaceMaterial() );
344
345 this->computeNd_matrixAt(ip->giveNaturalCoordinates(), N_d);
346 this->computeBd_matrixAt(ip, B_d);
347
348 // K_dd1 = ( 1 + l*f' ) * N^t*N;
349 // K_dd2 = ( l*l ) * B^t*B;
350 double dA = this->computeAreaAround(ip);
351
352 double fPrime = mat->giveDrivingForcePrime(ip);
353 double factorN = 1.0 + l* fPrime;
354 double factorB = l*l;
355
356 answer.plusProductUnsym(N_d, N_d, factorN * dA);
357 answer.plusProductUnsym(B_d, B_d, factorB * dA);
358 }
359
360 auto temp = answer;
361 answer.resize(4,4);
362 answer.zero();
363 answer.assemble(temp,{1,2},{1,2});
364 answer.assemble(temp,{3,4},{3,4});
365}
366
367
368void
369StructuralInterfaceElementPhF :: computeTraction(FloatArray &traction, IntegrationPoint *ip, FloatArray &jump, TimeStep *tStep)
370{
371 // Returns the traction in global coordinate system
372 FloatMatrix rotationMatGtoL, F;
373 this->computeTransformationMatrixAt(ip, rotationMatGtoL);
374 jump.rotatedWith(rotationMatGtoL, 'n'); // transform jump to local coord system
375
376 double damage = this->computeDamageAt(ip, VM_Total, tStep);
377
378 this->giveEngTraction(traction, ip, jump, damage, tStep);
379
380 traction.rotatedWith(rotationMatGtoL, 't'); // transform traction to global coord system
381}
382
383
384void
385StructuralInterfaceElementPhF :: updateYourself(TimeStep *tStep)
386{
387 Element :: updateYourself(tStep);
388
389 // record initial displacement if element not active
390 if ( activityTimeFunction && !isActivated(tStep) ) {
391 this->computeVectorOf(VM_Total, tStep, initialDisplacements);
392 }
393}
394
395
396void
397StructuralInterfaceElementPhF :: updateInternalState(TimeStep *tStep)
398{
399 // Updates the receiver at end of step.
400 FloatArray tractionG, jumpL;
401
402 // force updating strains & stresses
403 for ( auto &iRule: integrationRulesArray ) {
404 for ( auto &gp: *iRule ) {
405 this->computeSpatialJump(jumpL, gp, tStep);
406 this->computeTraction(tractionG, gp, jumpL, tStep);
407 }
408 }
409}
410
411
412} // end namespace oofem
413
#define N(a, b)
bool hasDofID(DofIDItem id) const
Definition dofmanager.C:174
int giveNumberOfDofs() const
Definition dofmanager.C:287
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
virtual FEInterpolation * giveInterpolation() const
Definition element.h:648
int activityTimeFunction
Element activity time function. If defined, nonzero value indicates active receiver,...
Definition element.h:163
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
virtual int giveNumberOfDofManagers() const
Definition element.h:695
DofManager * giveDofManager(int i) const
Definition element.C:553
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
void assemble(const FloatArray &fe, const IntArray &loc)
Definition floatarray.C:616
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Definition floatarray.C:288
double dotProduct(const FloatArray &x) const
Definition floatarray.C:524
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void rotatedWith(FloatMatrix &r, char mode)
Definition floatarray.C:814
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
void add(const FloatArray &src)
Definition floatarray.C:218
void subtract(const FloatArray &src)
Definition floatarray.C:320
void rotatedWith(const FloatMatrix &r, char mode='n')
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 beProductTOf(const FloatMatrix &a, const FloatMatrix &b)
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void zero()
Zeroes all coefficient of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
void followedBy(const IntArray &b, int allocChunk=0)
Definition intarray.C:94
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
void computeStiffnessMatrix_uu(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
virtual void computeCovarBaseVectorsAt(GaussPoint *gp, FloatMatrix &G)=0
virtual void getLocationArray_d(IntArray &answer)=0
virtual void computeNd_matrixAt(const FloatArray &lCoords, FloatMatrix &N)
void giveInternalForcesVector_d(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
virtual void getLocationArray_u(IntArray &answer)=0
virtual void computeTraction(FloatArray &traction, IntegrationPoint *ip, FloatArray &jump, TimeStep *tStep)
void giveInternalForcesVector_u(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
void computeStiffnessMatrix_dd(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
virtual void giveDofManDofIDMask_u(IntArray &answer)=0
virtual void giveEngTraction(FloatArray &answer, GaussPoint *gp, const FloatArray &jump, double damage, TimeStep *tStep)
virtual void giveDofManDofIDMask_d(IntArray &answer)=0
virtual double computeDamageAt(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
void computeDamageUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *stepN)
virtual void computeBd_matrixAt(GaussPoint *, FloatMatrix &)
FloatArray initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
virtual void computeTransformationMatrixAt(GaussPoint *gp, FloatMatrix &answer)=0
StructuralInterfaceCrossSection * giveInterfaceCrossSection()
virtual void computeSpatialJump(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
virtual double computeAreaAround(GaussPoint *gp)=0
virtual void giveStiffnessMatrix_Eng(FloatMatrix &answer, MatResponseMode rMode, IntegrationPoint *ip, TimeStep *tStep)
virtual void computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer)=0
virtual double giveDrivingForcePrime(GaussPoint *gp) const
virtual double giveDrivingForce(GaussPoint *gp) const
#define S(p)
Definition mdm.C:469
GaussPoint IntegrationPoint
Definition gausspoint.h:272

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