OOFEM 3.0
Loading...
Searching...
No Matches
transportelement.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
36#include "domain.h"
38#include "load.h"
39#include "bodyload.h"
40#include "boundaryload.h"
41#include "gausspoint.h"
43#include "intarray.h"
44#include "floatarray.h"
45#include "floatmatrix.h"
46#include "verbose.h"
47#include "mathfem.h"
48#include "crosssection.h"
50#include "feinterpol.h"
51#include "feinterpol1d.h"
52#include "feinterpol2d.h"
53#include "feinterpol3d.h"
54#include "dof.h"
56#include "function.h"
57#include "parametermanager.h"
58#include "paramkey.h"
59
60#ifdef __CEMHYD_MODULE
62#endif
63
64
65#ifdef __OOFEG
66 #include "oofeggraphiccontext.h"
67#endif
68
69namespace oofem {
70
71const double TransportElement :: stefanBoltzmann = 5.67e-8; //W/m2/K4
73
74
75TransportElement :: TransportElement(int n, Domain *aDomain, ElementMode em) :
76 Element(n, aDomain), emode( em )
77{
78 this->vofFunction = 0;
79}
80
81
82void
83TransportElement :: initializeFrom(InputRecord &ir, int priority)
84{
85 Element::initializeFrom(ir, priority);
86 ParameterManager &ppm = this->giveDomain()->elementPPM;
88}
89
90void
91TransportElement :: postInitialize()
92{
94}
95
96void
97TransportElement :: giveDofManDofIDMask(int inode, IntArray &answer) const
98{
99 if ( emode == HeatTransferEM ) {
100 answer = {T_f};
101 } else if ( emode == HeatMass1TransferEM ) {
102 answer = {T_f, C_1};
103 } else if ( emode == Mass1TransferEM ) {
104 answer = {C_1};
105 } else {
106 OOFEM_ERROR("Unknown ElementMode");
107 }
108}
109
110
111void
112TransportElement :: giveCharacteristicMatrix(FloatMatrix &answer,
113 CharType mtrx, TimeStep *tStep)
114//
115// returns characteristics matrix of receiver according to mtrx
116//
117{
118 if ( mtrx == TangentStiffnessMatrix ) {
119 FloatMatrix tmp;
120 this->computeConductivityMatrix(answer, Conductivity, tStep);
121 this->computeBCMtrxAt(tmp, tStep, VM_TotalIntrinsic);
122 answer.add(tmp);
123 } else if ( mtrx == ConductivityMatrix ) {
124 this->computeConductivityMatrix(answer, Conductivity, tStep);
125 } else if ( mtrx == CapacityMatrix || mtrx == MassMatrix ) {
126 this->computeCapacityMatrix(answer, tStep);
127 } else if ( mtrx == LumpedMassMatrix ) {
128 this->computeCapacityMatrix(answer, tStep);
129 for ( int i = 1; i <= answer.giveNumberOfRows(); i++ ) {
130 double s = 0.0;
131 for ( int j = 1; j <= answer.giveNumberOfColumns(); j++ ) {
132 s += answer.at(i, j);
133 answer.at(i, j) = 0.0;
134 }
135 answer.at(i, i) = s;
136 }
137 } else {
138 OOFEM_ERROR("Unknown Type of characteristic mtrx (%s)", __CharTypeToString(mtrx));
139 }
140}
141
142
143void
144TransportElement :: giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode,
145 TimeStep *tStep)
146//
147// returns characteristics vector of receiver according to requested type
148//
149{
150 if ( mtrx == InternalForcesVector ) {
151 this->computeInternalForcesVector(answer, tStep);
152 } else if ( mtrx == ExternalForcesVector ) {
153 //this->computeExternalForcesVector(answer, tStep, mode); // bp: assembled by emodel->updateComponent
154 answer.resize(0);
155 } else if ( mtrx == InertiaForcesVector ) {
156 this->computeInertiaForcesVector(answer, tStep);
157 } else if ( mtrx == LumpedMassMatrix ) {
158 this->computeLumpedCapacityVector(answer, tStep);
159 } else {
160 OOFEM_ERROR( "Unknown Type of characteristic mtrx (%s)",
161 __CharTypeToString(mtrx) );
162 }
163}
164
165
166int
167TransportElement :: checkConsistency()
168//
169// check internal consistency
170// mainly tests, whether material and crossSection data
171// are safe for conversion to "Structural" versions
172//
173{
174 int result = 1;
175 if ( this->material > 0 && !dynamic_cast< TransportMaterial* >( this->giveMaterial() ) ) {
176 OOFEM_WARNING("cross section without support for transport problems");
177 result = 0;
178 }
179
180#if 0
181 if ( !this->giveCrossSection()->testCrossSectionExtension(CS_TransportCapability) ) {
182 OOFEM_WARNING("cross-section without support for transport problems", 1);
183 result = 0;
184 }
185
186#endif
187 return result;
188}
189
190void
191TransportElement :: computeCapacityMatrix(FloatMatrix &answer, TimeStep *tStep)
192{
193 answer.resize( this->computeNumberOfDofs(), this->computeNumberOfDofs() );
194 answer.zero();
195
196 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
197 this->computeCapacitySubMatrix(answer, Capacity, 0, tStep);
198 } else if ( emode == HeatMass1TransferEM ) {
199 FloatMatrix subAnswer;
200 MatResponseMode rmode [ 2 ] = {
201 Capacity_hh, Capacity_ww
202 };
203 //double coeff = 1.0; //this->giveMaterial()->give('d');
204
205 for ( int i = 1; i <= 2; i++ ) {
206 this->computeCapacitySubMatrix(subAnswer, rmode [ i - 1 ], 0, tStep);
207 this->assembleLocalContribution(answer, subAnswer, 2, i, i);
208 }
209 } else {
210 OOFEM_ERROR("Unknown ElementMode");
211 }
212}
213
214void
215TransportElement :: computeConductivityMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
216{
217 answer.resize( this->computeNumberOfDofs(), this->computeNumberOfDofs() );
218 answer.zero();
219 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
220 this->computeConductivitySubMatrix(answer, 0, Conductivity_hh, tStep);
221 } else if ( emode == HeatMass1TransferEM ) {
222 FloatMatrix subAnswer;
223 MatResponseMode rmode [ 2 ] [ 2 ] = { { Conductivity_hh, Conductivity_hw }, { Conductivity_wh, Conductivity_ww } };
224
225 for ( int i = 1; i <= 2; i++ ) {
226 for ( int j = 1; j <= 2; j++ ) {
227 this->computeConductivitySubMatrix(subAnswer, 0, rmode [ i - 1 ] [ j - 1 ], tStep);
228 this->assembleLocalContribution(answer, subAnswer, 2, i, j);
229 }
230 }
231 } else {
232 OOFEM_ERROR("Unknown ElementMode");
233 }
234}
235
236
237void
238TransportElement :: computeNAt(FloatArray &answer, const FloatArray &lcoord)
239{
240 this->giveInterpolation()->evalN( answer, lcoord, FEIElementGeometryWrapper(this) );
241}
242
243
244void
245TransportElement :: computeNmatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
246{
247 FloatArray n;
248 this->computeNAt(n, lcoords);
249 if ( this->emode == HeatTransferEM || this->emode == Mass1TransferEM ) {
250 answer.beNMatrixOf(n, 1);
251 } else {
252 answer.beNMatrixOf(n, 2);
253 }
254}
255
256
257void
258TransportElement :: computeBmatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
259{
260 FloatMatrix dnx;
262 this->giveInterpolation()->evaldNdx( dnx, lcoords, FEIElementGeometryWrapper(this) );
263 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
264 answer.beTranspositionOf(dnx);
265 } else if ( this->emode == HeatMass1TransferEM ) {
266 int nodes = dnx.giveNumberOfRows();
267 int nsd = dnx.giveNumberOfColumns();
268 answer.resize(nsd*2, nodes*2);
269 answer.zero();
270 for (int i = 0; i < nodes; ++i) {
271 for (int j = 0; j < nsd; ++j) {
272 answer(j, i*2) = dnx(i, j);
273 answer(j+nsd, i*2+1) = dnx(i, j);
274 }
275 }
276 }
277}
278
279
280void
281TransportElement :: computeGradientMatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
282{
283 FloatMatrix dnx;
284 this->giveInterpolation()->evaldNdx( dnx, lcoords, FEIElementGeometryWrapper(this) );
286 answer.beTranspositionOf(dnx);
287}
288
289
290void
291TransportElement :: computeEgdeNAt(FloatArray &answer, int iedge, const FloatArray &lcoords)
292{
293 FEInterpolation *interp = this->giveInterpolation();
294 if ( dynamic_cast< FEInterpolation2d * >(interp) ) {
295 dynamic_cast< FEInterpolation2d * >(interp)->edgeEvalN( answer, iedge, lcoords, FEIElementGeometryWrapper(this) );
296 } else if ( dynamic_cast< FEInterpolation3d * >(interp) ) {
297 dynamic_cast< FEInterpolation3d * >(interp)->edgeEvalN( answer, iedge, lcoords, FEIElementGeometryWrapper(this) );
298 }
299}
300
301
302void
303TransportElement :: giveEdgeDofMapping(IntArray &answer, int iEdge)
304{
305 FEInterpolation *interp = this->giveInterpolation();
306 if ( dynamic_cast< FEInterpolation2d * >(interp) ) {
307 answer = dynamic_cast< FEInterpolation2d * >(interp)->computeLocalEdgeMapping(iEdge);
308 } else if ( dynamic_cast< FEInterpolation3d * >(interp) ) {
309 answer = dynamic_cast< FEInterpolation3d * >(interp)->computeLocalEdgeMapping(iEdge);
310 }
311}
312
313void
314TransportElement :: computeSurfaceNAt(FloatArray &answer, int iSurf, const FloatArray &lcoord)
315{
316 FEInterpolation *interp = this->giveInterpolation();
317 if ( dynamic_cast< FEInterpolation3d * >(interp) ) {
318 dynamic_cast< FEInterpolation3d * >(interp)->surfaceEvalN( answer, iSurf, lcoord, FEIElementGeometryWrapper(this) );
319 }
320}
321
322void
323TransportElement :: giveSurfaceDofMapping(IntArray &answer, int iSurf)
324{
325 FEInterpolation *interp = this->giveInterpolation();
326 if ( dynamic_cast< FEInterpolation3d * >(interp) ) {
327 answer = dynamic_cast< FEInterpolation3d * >(interp)->computeLocalSurfaceMapping(iSurf);
328 }
329}
330
331int
332TransportElement :: giveApproxOrder(int unknownIndx)
333{
334 return this->giveInterpolation()->giveInterpolationOrder();
335}
336
337void
338TransportElement :: computeCapacitySubMatrix(FloatMatrix &answer, MatResponseMode rmode, int iri, TimeStep *tStep)
339{
340 FloatArray n;
341 TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
342
343 answer.clear();
344 for ( auto &gp: *integrationRulesArray [ iri ] ) {
345 this->computeNAt( n, gp->giveNaturalCoordinates() );
346 // ask for capacity coefficient. In basic units [J/K/m3]
347 double c = mat->giveCharacteristicValue(rmode, gp, tStep);
348 double dV = this->computeVolumeAround(gp);
349 answer.plusDyadSymmUpper(n, dV * c);
350 }
351
352 answer.symmetrized();
353}
354
355void
356TransportElement :: computeConductivitySubMatrix(FloatMatrix &answer, int iri, MatResponseMode rmode, TimeStep *tStep)
357{
358 FloatMatrix b, d, db;
359
360 answer.resize( this->giveNumberOfDofManagers(), this->giveNumberOfDofManagers() );
361 answer.zero();
362 for ( auto &gp: *integrationRulesArray [ iri ] ) {
363 this->computeConstitutiveMatrixAt(d, rmode, gp, tStep);
364 this->computeGradientMatrixAt(b, gp->giveNaturalCoordinates());
365 double dV = this->computeVolumeAround(gp);
366
367 db.beProductOf(d, b);
368 answer.plusProductSymmUpper(b, db, dV);
369 }
370
371 answer.symmetrized();
372}
373
374void
375TransportElement :: computeInternalSourceRhsSubVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode, int indx)
376{
377 // Computes numerically the generator Rhs vector of the receiver due to the generator
378 // at tStep.
379 // // load is first transformed to local cs.
380 // // load vector is then transformed to coordinate system in each node.
381 // // (should be global coordinate system, but there may be defined
382 // // different coordinate system in each node)
383 TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
384
385 FloatArray val, n;
386 answer.clear();
387
388 // add internal source produced by material (if any)
389 if ( mat->hasInternalSource() ) {
390 for ( auto &gp: *this->giveDefaultIntegrationRulePtr() ) {
391 this->computeNAt( n, gp->giveNaturalCoordinates() );
392 double dV = this->computeVolumeAround(gp);
393 mat->computeInternalSourceVector(val, gp, tStep, mode);
394
395 answer.add(val.at(indx) * dV, n);
396 }
397 }
398}
399
400
401void
402TransportElement :: computeInternalSourceRhsVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
403{
404 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
405 this->computeInternalSourceRhsSubVectorAt(answer, tStep, mode, 1);
406 } else if ( emode == HeatMass1TransferEM ) {
407 FloatArray subAnswer;
408
409 for ( int i = 1; i <= 2; i++ ) {
410 this->computeInternalSourceRhsSubVectorAt(subAnswer, tStep, mode, i);
411 if ( subAnswer.isNotEmpty() ) {
412 if ( answer.isEmpty() ) {
413 answer.resize( 2 * subAnswer.giveSize() );
414 answer.zero();
415 }
416
417 this->assembleLocalContribution(answer, subAnswer, 2, i);
418 }
419 }
420 } else {
421 OOFEM_ERROR("Unknown ElementMode");
422 }
423}
424
425
426void
427TransportElement :: computeIntSourceLHSMatrix(FloatMatrix &answer, TimeStep *tStep)
428{
429 TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
430 if ( mat->hasInternalSource() ) {
431 answer.resize( this->computeNumberOfDofs(), this->computeNumberOfDofs() );
432 answer.zero();
433
434 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
435 this->computeIntSourceLHSSubMatrix(answer, IntSource, 0, tStep);
436 } else if ( emode == HeatMass1TransferEM ) {
437 FloatMatrix subAnswer;
438 MatResponseMode rmode [ 2 ] = {
439 IntSource_hh, IntSource_ww
440 };
441 //double coeff = 1.0; //this->giveMaterial()->give('d');
442
443 for ( int i = 1; i <= 2; i++ ) {
444 this->computeIntSourceLHSSubMatrix(subAnswer, rmode [ i - 1 ], 0, tStep);
445 this->assembleLocalContribution(answer, subAnswer, 2, i, i);
446 }
447 } else {
448 OOFEM_ERROR("Unknown ElementMode");
449 }
450 } else {
451 answer.clear();
452 }
453}
454
455void
456TransportElement :: computeIntSourceLHSSubMatrix(FloatMatrix &answer, MatResponseMode rmode, int iri, TimeStep *tStep)
457// computes LHS matrix due to material internal source (dHeat/dTemperature, dWaterSource/dw)
458// IntSource - Heat transfer
459// IntSource_hh - HeMo heat source
460// IntSource_ww - HeMo water content source
461// hw, wh is not used now
462{
463 FloatArray n;
464 TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
465
466 answer.clear();
467 for ( auto &gp: *integrationRulesArray [ iri ] ) {
468 this->computeNAt( n, gp->giveNaturalCoordinates() );
469 // ask for coefficient from material
470 double c = mat->giveCharacteristicValue(rmode, gp, tStep);
471 double dV = this->computeVolumeAround(gp);
472 answer.plusDyadSymmUpper(n, dV * c);
473 }
474
475 answer.symmetrized();
476}
477
478
479void
480TransportElement :: computeConstitutiveMatrixAt(FloatMatrix &answer,
481 MatResponseMode rMode, GaussPoint *gp,
482 TimeStep *tStep)
483{
484 static_cast< TransportMaterial * >( this->giveMaterial() )->giveCharacteristicMatrix(answer, rMode, gp, tStep);
485}
486
487void
488TransportElement :: computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
489{
490 FloatArray unknowns;
491 this->computeVectorOf(VM_TotalIntrinsic, tStep, unknowns);
492
493 TransportMaterial *mat = static_cast< TransportMaterial* >( this->giveMaterial() );
494 FloatArray flux, grad, field;
495 FloatMatrix B, N;
496
497 answer.clear();
498 for ( auto &gp: *integrationRulesArray [ 0 ] ) {
499 const FloatArray &lcoords = gp->giveNaturalCoordinates();
500
501 this->computeNmatrixAt(N, lcoords);
502 this->computeBmatrixAt(B, lcoords);
503
504 field.beProductOf(N, unknowns);
505 grad.beProductOf(B, unknowns);
506
507 mat->giveFluxVector(flux, gp, grad, field, tStep);
508
509 double dV = this->computeVolumeAround(gp);
510 answer.plusProduct(B, flux, -dV);
511
513 if ( mat->hasInternalSource() ) {
514 // add internal source produced by material (if any)
515 FloatArray val;
516 mat->computeInternalSourceVector(val, gp, tStep, VM_TotalIntrinsic);
517 answer.plusProduct(N, val, -dV);
518 }
519 }
520 /*
522 // Neumann b.c.s
523 FloatMatrix bc_tangent;
524 this->computeBCMtrxAt(bc_tangent, tStep, VM_TotalIntrinsic);
525 if ( bc_tangent.isNotEmpty() ) {
526 answer.plusProduct(bc_tangent, unknowns, 1.0);
527 }
528 */
529}
530
531
532void
533TransportElement :: computeInertiaForcesVector(FloatArray &answer, TimeStep *tStep)
534{
535 FloatMatrix cap;
536 FloatArray vel;
537 this->computeVectorOf(VM_Velocity, tStep, vel);
538 this->computeCapacityMatrix(cap, tStep);
539 answer.beProductOf(cap, vel);
540}
541
542
543void
544TransportElement :: computeLumpedCapacityVector(FloatArray &answer, TimeStep *tStep)
545{
546 FloatMatrix cap;
547 this->computeCapacityMatrix(cap, tStep);
548
549 int size = cap.giveNumberOfRows();
550 answer.resize(size);
551 answer.zero();
552 for ( int i = 1; i <= size; i++ ) {
553 for ( int j = 1; j <= size; j++ ) {
554 answer.at(i) += cap.at(i, j);
555 }
556 }
557}
558
559
560void
561TransportElement :: computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
562{
563 answer.clear();
564
566#if 1
567 if ( type != ExternalForcesVector ) {
568 return;
569 }
570#else
571 if ( !( load->giveType() == TransmissionBC && type == ExternalForcesVector ) &&
572 !( load->giveType() == ConvectionBC && type == InternalForcesVector ) ) {
573 return;
574 }
575#endif
576
577 FloatArray gcoords, val, n, unknowns;
579 IntArray dofid;
580 int unknownsPerNode = this->emode == HeatMass1TransferEM ? 2 : 1;
581
582 this->giveElementDofIDMask(dofid);
583
585 FEInterpolation *interp = this->giveInterpolation();
586 std :: unique_ptr< IntegrationRule > iRule( interp->giveIntegrationRule( load->giveApproxOrder() + 1 + interp->giveInterpolationOrder(), this->giveGeometryType()) );
587 iRule->setElement(this);
588
589 if ( load->giveType() == ConvectionBC || load->giveType() == RadiationBC ) {
590 this->computeVectorOf(dofid, VM_TotalIntrinsic, tStep, unknowns);
591 }
592
593 for ( auto &gp: *iRule ) {
594 const FloatArray &lcoords = gp->giveNaturalCoordinates();
595
596 interp->evalN( n, lcoords, FEIElementGeometryWrapper(this) );
597 N.beNMatrixOf(n, unknownsPerNode);
598
599 double dV = this->computeVolumeAround(gp);
600 load->computeValueAt(val, tStep, gp, mode);
601 answer.plusProduct(N, val, dV);
602 }
603}
604
605
606void
607TransportElement :: computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
608{
609 answer.clear();
610
611 if ( !( load->giveType() == TransmissionBC && type == ExternalForcesVector ) &&
612 !( load->giveType() == ConvectionBC && type == InternalForcesVector ) &&
613 !( load->giveType() == RadiationBC && type == InternalForcesVector ) ) {
614 return;
615 }
616
617 FloatArray gcoords, val, n, unknowns, field;
619 IntArray dofid;
620 int unknownsPerNode = this->emode == HeatMass1TransferEM ? 2 : 1;
621
622 this->giveElementDofIDMask(dofid);
623
625 FEInterpolation *interp = this->giveInterpolation();
626 std :: unique_ptr< IntegrationRule > iRule( this->giveBoundarySurfaceIntegrationRule(load->giveApproxOrder() + interp->giveInterpolationOrder(), boundary) );
627
628 if ( load->giveType() == ConvectionBC || load->giveType() == RadiationBC ) {
629 const auto &bNodes = interp->boundaryGiveNodes(boundary, this->giveGeometryType());
630 this->computeBoundaryVectorOf(bNodes, dofid, VM_TotalIntrinsic, tStep, unknowns);
631 }
632
633 for ( auto &gp: *iRule ) {
634 const FloatArray &lcoords = gp->giveNaturalCoordinates();
635
636 //this->computeEgdeNAt(n, iEdge, lcoords);
637 interp->boundaryEvalN( n, boundary, lcoords, FEIElementGeometryWrapper(this) );
638 N.beNMatrixOf(n, unknownsPerNode);
639
640 interp->boundaryLocal2Global( gcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
641 double detJ = interp->boundaryGiveTransformationJacobian( boundary, lcoords, FEIElementGeometryWrapper(this) );
642 double dA = this->giveThicknessAt(gcoords) * gp->giveWeight() * detJ;
643
644 //Check external ambient temperature field first
645 FieldPtr tf;
646 if ((tf = domain->giveEngngModel()->giveContext()->giveFieldManager()->giveField(FT_TemperatureAmbient))){
647 tf->evaluateAt(val, gcoords, VM_TotalIntrinsic, tStep);
648 } else if ( load->giveFormulationType() == Load :: FT_Entity ) {
649 load->computeValueAt(val, tStep, lcoords, mode);
650 } else {
651 load->computeValueAt(val, tStep, gcoords, mode);
652 }
653
654 if ( load->giveType() == TransmissionBC ) {
655 val.negated();
656 } else if ( load->giveType() == ConvectionBC ) {
657 field.beProductOf(N, unknowns);
658 val.subtract(field);
659 double a;
660 if ( load->propertyMultExpr.isDefined() ) {
661 if ( unknownsPerNode != 1 ) {
662 OOFEM_ERROR("Not implemented for >=2 coupled fields");
663 }
664 a = load->giveProperty('a', tStep, { { "u", field.at(1) } });
665 } else {
666 a = load->giveProperty('a', tStep);
667 }
668 val.times( -1.0 * a );
669 } else if ( load->giveType() == RadiationBC ) {
670 //actual Temperature in C in field
671 field.beProductOf(N, unknowns); field.add(load->giveTemperOffset()); field.power(4);
672 val.add(load->giveTemperOffset()); val.power(4); val.subtract(field);
673 val.times( -1.0 * load->giveProperty('e', tStep) * stefanBoltzmann );
674 }
675
676 answer.plusProduct(N, val, dA);
677 }
678}
679
680//Contribution to conductivity matrix from convection and radiation
681void
682TransportElement :: computeTangentFromSurfaceLoad(FloatMatrix &answer, BoundaryLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep)
683{
684 answer.clear();
685
686 if ( load->giveType() != ConvectionBC ) {
687 return;
688 }
689
690 FloatArray gcoords, n, unknowns;
692 int unknownsPerNode = this->emode == HeatMass1TransferEM ? 2 : 1;
693
695 FEInterpolation *interp = this->giveInterpolation();
696 std :: unique_ptr< IntegrationRule > iRule( this->giveBoundarySurfaceIntegrationRule(load->giveApproxOrder() + interp->giveInterpolationOrder(), boundary) );
697
698 if ( load->propertyMultExpr.isDefined() ) {
699 if ( unknownsPerNode != 1 ) {
700 OOFEM_ERROR("Load property multexpr not implemented for coupled fields");
701 }
702 IntArray dofid;
703 const auto &bNodes = interp->boundaryGiveNodes(boundary, this->giveGeometryType());
704 this->giveElementDofIDMask(dofid);
705 this->computeBoundaryVectorOf(bNodes, dofid, VM_TotalIntrinsic, tStep, unknowns);
706 }
707
708 for ( auto &gp : *iRule ) {
709 const FloatArray &lcoords = gp->giveNaturalCoordinates();
710
711 interp->boundaryEvalN( n, boundary, lcoords, FEIElementGeometryWrapper(this) );
712 interp->boundaryLocal2Global( gcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
713 double detJ = interp->boundaryGiveTransformationJacobian( boundary, lcoords, FEIElementGeometryWrapper(this) );
714 double dA = this->giveThicknessAt(gcoords) * gp->giveWeight() * detJ;
715 N.beNMatrixOf(n, unknownsPerNode);
716
717 if ( load->giveType() == ConvectionBC ){
718 double a;
719 if ( load->propertyMultExpr.isDefined() ) {
720 double value = n.dotProduct(unknowns);
721 a = load->giveProperty('a', tStep, { { "u", value } });
722 } else {
723 a = load->giveProperty('a', tStep);
724 }
725 answer.plusProductSymmUpper(N, N, a * dA);
726 }
727 }
728 answer.symmetrized();
729}
730
731
732void
733TransportElement :: computeTangentFromEdgeLoad(FloatMatrix &answer, BoundaryLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep)
734{
735 answer.clear();
736
737 if ( load->giveType() != ConvectionBC ) {
738 return;
739 }
740
741 FloatArray gcoords, n, unknowns;
743 int unknownsPerNode = this->emode == HeatMass1TransferEM ? 2 : 1;
744
746 FEInterpolation *interp = this->giveInterpolation();
747 std :: unique_ptr< IntegrationRule > iRule( this->giveBoundaryEdgeIntegrationRule(load->giveApproxOrder() + interp->giveInterpolationOrder(), boundary) );
748
749 if ( load->propertyMultExpr.isDefined() ) {
750 IntArray dofid;
751 this->giveElementDofIDMask(dofid);
752 const auto &bNodes = interp->boundaryEdgeGiveNodes(boundary, this->giveGeometryType());
753 this->giveElementDofIDMask(dofid);
754 this->computeBoundaryVectorOf(bNodes, dofid, VM_TotalIntrinsic, tStep, unknowns);
755 }
756
757 for ( auto &gp : *iRule ) {
758 const FloatArray &lcoords = gp->giveNaturalCoordinates();
759
760 interp->boundaryEvalN( n, boundary, lcoords, FEIElementGeometryWrapper(this) );
761 interp->boundaryLocal2Global( gcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
762 double detJ = interp->boundaryGiveTransformationJacobian( boundary, lcoords, FEIElementGeometryWrapper(this) );
763 double dA = this->giveThicknessAt(gcoords) * gp->giveWeight() * detJ;
764 N.beNMatrixOf(n, unknownsPerNode);
765
766 if ( load->giveType() == ConvectionBC ){
767 double a;
768 if ( load->propertyMultExpr.isDefined() ) {
769 if ( unknownsPerNode != 1 ) {
770 OOFEM_ERROR("Not implemented for >=2 coupled fields");
771 }
772 FloatArray field, val;
773 field.beProductOf(N, unknowns);
774 val.subtract(field);
775 a = load->giveProperty('a', tStep, { { "u", field.at(1) } });
776 } else {
777 a = load->giveProperty('a', tStep);
778 }
779 answer.plusProductSymmUpper(N, N, a * dA);
780 }
781 }
782 answer.symmetrized();
783}
784
785
786void
787TransportElement :: computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
788{
790 answer.clear();
791
792 if ( !( load->giveType() == TransmissionBC && type == ExternalForcesVector ) &&
793 !( load->giveType() == ConvectionBC && type == InternalForcesVector ) &&
794 !( load->giveType() == RadiationBC && type == InternalForcesVector ) ) {
795 return;
796 }
797
798 if ( !load->isImposed(tStep) ) {
799 return;
800 }
801
802 FloatArray gcoords, val, n, unknowns, field;
804 int unknownsPerNode = this->emode == HeatMass1TransferEM ? 2 : 1;
805
807 FEInterpolation *interp = this->giveInterpolation();
808 std :: unique_ptr< IntegrationRule > iRule( this->giveBoundaryEdgeIntegrationRule(load->giveApproxOrder() + interp->giveInterpolationOrder(), boundary) );
809
810 // get solution vector for InternalForcesVector (Convention or radiation)
811 if ( load->giveType() == ConvectionBC || load->giveType() == RadiationBC ) {
812 IntArray dofid;
813 this->giveElementDofIDMask(dofid);
814 const auto &bNodes = interp->boundaryEdgeGiveNodes(boundary, this->giveGeometryType());
815 this->computeBoundaryVectorOf(bNodes, dofid, VM_TotalIntrinsic, tStep, unknowns);
816 }
817
818 for ( auto &gp: *iRule ) {
819 const FloatArray &lcoords = gp->giveNaturalCoordinates();
820 //this->computeEgdeNAt(n, boundary, lcoords);
821 interp->boundaryEdgeEvalN( n, boundary, lcoords, FEIElementGeometryWrapper(this) );
822 N.beNMatrixOf(n, unknownsPerNode);
823
824
825 //double detJ = interp->boundaryEdgeGiveTransformationJacobian( boundary, lcoords, FEIElementGeometryWrapper(this) );
826 //double dL = this->giveThicknessAt(gcoords) * gp->giveWeight() * detJ;
827 double dV = this->computeEdgeVolumeAround(gp, boundary);
828
829 FieldPtr tf = domain->giveEngngModel()->giveContext()->giveFieldManager()->giveField(FT_TemperatureAmbient);
830 if ( tf ) {
831 interp->boundaryLocal2Global( gcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
832 tf->evaluateAt(val, gcoords, VM_TotalIntrinsic, tStep);
833 } else if ( load->giveFormulationType() == Load :: FT_Entity ) {
834 load->computeValueAt(val, tStep, lcoords, mode);
835 } else {
836 interp->boundaryLocal2Global( gcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
837 load->computeValueAt(val, tStep, gcoords, mode);
838 }
839
840 if ( load->giveType() == TransmissionBC ) {
841 val.negated();
842 } else if ( load->giveType() == ConvectionBC ) {
843 field.beProductOf(N, unknowns);
844 val.subtract(field);
845 double a;
846
847 if ( load->propertyMultExpr.isDefined() ) {
848 if ( unknownsPerNode != 1 ) {
849 OOFEM_ERROR("Not implemented for >=2 coupled fields");
850 }
851 a = load->giveProperty('a', tStep, { { "u", field.at(1) } });
852 } else {
853 a = load->giveProperty('a', tStep);
854 }
855 val.times( -1.0 * a );
856 } else if ( load->giveType() == RadiationBC ) {
857 //actual Temperature in C in field
858 field.beProductOf(N, unknowns); field.add(load->giveTemperOffset()); field.power(4);
859 val.add(load->giveTemperOffset()); val.power(4); val.subtract(field);
860 val.times( -1.0 * load->giveProperty('e', tStep) * stefanBoltzmann );
861 }
862 answer.plusProduct(N, val, dV);
863 }
864}
865
866void
867TransportElement :: computeExternalForcesVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
868{
869 this->computeBCVectorAt(answer, tStep, mode);
870}
871
872void
873TransportElement :: computeBCVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
874{
875 answer.resize( this->computeNumberOfDofs() );
876 answer.zero();
877
878 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
879 this->computeBCSubVectorAt(answer, tStep, mode, 1);
880 } else if ( emode == HeatMass1TransferEM ) {
881 FloatArray subAnswer;
882
883 for ( int i = 1; i <= 2; i++ ) {
884 this->computeBCSubVectorAt(subAnswer, tStep, mode, i);
885 this->assembleLocalContribution(answer, subAnswer, 2, i);
886 }
887 } else {
888 OOFEM_ERROR("Unknown ElementMode");
889 }
890}
891
892
893void
894TransportElement :: computeBCMtrxAt(FloatMatrix &answer, TimeStep *tStep, ValueModeType mode)
895{
896 int ndofs = this->computeNumberOfDofs();
897 answer.resize(ndofs, ndofs);
898 answer.zero();
899
900 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
901 this->computeBCSubMtrxAt(answer, tStep, mode, 1);
902 } else if ( emode == HeatMass1TransferEM ) {
903 FloatMatrix subAnswer;
904
905 for ( int i = 1; i <= 2; i++ ) {
906 this->computeBCSubMtrxAt(subAnswer, tStep, mode, i);
907 if ( subAnswer.isNotEmpty() ) {
908 this->assembleLocalContribution(answer, subAnswer, 2, i, i);
909 }
910 }
911 } else {
912 OOFEM_ERROR("Unknown ElementMode");
913 }
914}
915
916
917void
918TransportElement :: computeBCSubVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode, int indx)
919{
920 FloatArray vec;
921
922 answer.resize( this->giveNumberOfDofManagers() );
923 answer.zero();
924
925 // loop over boundary load array
926 int nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
927 for ( int i = 1; i <= nLoads; i++ ) {
928 int n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
929 int id = boundaryLoadArray.at(i * 2);
930 Load *load = domain->giveLoad(n);
931 if ( load->isDofExcluded(indx) || !load->isImposed(tStep) ) {
932 continue;
933 }
934
935 bcGeomType ltype = load->giveBCGeoType();
936 if ( ltype == EdgeLoadBGT ) {
937 this->computeEdgeBCSubVectorAt(vec, load, id, tStep, mode, indx);
938 } else if ( ltype == SurfaceLoadBGT ) {
939 this->computeSurfaceBCSubVectorAt(vec, load, id, tStep, mode, indx);
940 } else {
941 OOFEM_ERROR("unsupported bc type encountered");
942 }
943
944 answer.add(vec);
945 }
946
947 for ( int i = 1; i <= bodyLoadArray.giveSize(); ++i ) {
948 Load *load = domain->giveLoad(bodyLoadArray.at(i));
949 if ( !load->isImposed(tStep) ) {
950 continue;
951 }
952 this->computeBodyBCSubVectorAt(answer, load, tStep, mode, indx);
953 }
954}
955
956
957void
958TransportElement :: computeBodyBCSubVectorAt(FloatArray &answer, Load *load,
959 TimeStep *tStep, ValueModeType mode, int indx)
960{
961 FloatArray val, globalIPcoords, n;
962 answer.resize( this->giveNumberOfDofManagers() );
963
964 std :: unique_ptr< IntegrationRule > iRule( this->giveInterpolation()->giveIntegrationRule(load->giveApproxOrder(), this->giveGeometryType()) );
965 for ( GaussPoint *gp : *iRule ) {
966 double dV = this->computeVolumeAround(gp);
967 this->computeNAt( n, gp->giveNaturalCoordinates() );
968 this->computeGlobalCoordinates( globalIPcoords, gp->giveNaturalCoordinates() );
969 load->computeValueAt(val, tStep, globalIPcoords, mode);
970 answer.add(val.at(indx) * dV, n);
971 }
972}
973
974
975void
976TransportElement :: computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge,
977 TimeStep *tStep, ValueModeType mode, int indx)
978{
979 answer.resize( this->giveNumberOfDofManagers() );
980 answer.zero();
981
982 if ( ( load->giveType() == TransmissionBC ) || ( load->giveType() == ConvectionBC ) || load->giveType() == RadiationBC ) {
983 BoundaryLoad *edgeLoad = static_cast< BoundaryLoad * >(load);
984
985 int approxOrder = edgeLoad->giveApproxOrder() + this->giveApproxOrder(indx);
986 int numberOfEdgeIPs = ( int ) ceil( ( approxOrder + 1.0 ) / 2. );
987 GaussIntegrationRule iRule(1, this, 1, 1);
988 iRule.SetUpPointsOnLine(numberOfEdgeIPs, _Unknown);
989 FloatArray reducedAnswer, val, n;
990 IntArray mask;
991 double dV, coeff = 1.0;
992
993 for ( GaussPoint *gp: iRule ) {
994 if ( load->giveType() == TransmissionBC ) {
995 coeff = -1.0;
996 } else if ( load->giveType() == ConvectionBC ) {
997 if( edgeLoad->propertyMultExpr.isDefined()) {//dependence on state variable
1000 OOFEM_ERROR("Not implemented for >=2 coupled fields");
1001 }
1002 FloatArray unknowns;
1003 IntArray dofid;
1004 this->computeEgdeNAt( n, iEdge, gp->giveNaturalCoordinates() );
1005 this->giveElementDofIDMask(dofid);
1006 this->giveEdgeDofMapping(mask, iEdge);
1007 this->computeBoundaryVectorOf(mask, dofid, VM_TotalIntrinsic, tStep, unknowns);
1008 double value = n.dotProduct(unknowns);
1009 coeff = edgeLoad->giveProperty('a', tStep, { { "u", value } });
1010 } else {
1011 coeff = edgeLoad->giveProperty('a', tStep);
1012 }
1013 } else if ( load->giveType() == RadiationBC ){
1014 OOFEM_ERROR("Not implemented, use TransientTransport solver");
1015 //coeff = getRadiativeHeatTranferCoef(edgeLoad, tStep);
1016 }
1017
1018 const FloatArray &lcoords = gp->giveNaturalCoordinates();
1019 this->computeEgdeNAt(n, iEdge, lcoords);
1020 dV = this->computeEdgeVolumeAround(gp, iEdge);
1021
1022 FieldPtr tf;
1023 FloatArray gcoords;
1024 if ((tf = domain->giveEngngModel()->giveContext()->giveFieldManager()->giveField(FT_TemperatureAmbient))){
1025 //this->computeEdgeIpGlobalCoords(gcoords, lcoords, iEdge);
1026 this->giveInterpolation()->boundaryEdgeLocal2Global( gcoords, iEdge, lcoords, FEIElementGeometryWrapper(this) );
1027
1028 tf->evaluateAt(val, gcoords, VM_TotalIntrinsic, tStep);
1029 } else if ( edgeLoad->giveFormulationType() == Load :: FT_Entity ) {
1030 edgeLoad->computeValueAt(val, tStep, lcoords, mode);
1031 } else {
1032 //this->computeEdgeIpGlobalCoords(gcoords, lcoords, iEdge);
1033 this->giveInterpolation()->boundaryEdgeLocal2Global( gcoords, iEdge, lcoords, FEIElementGeometryWrapper(this) );
1034 edgeLoad->computeValueAt(val, tStep, gcoords, mode);
1035 }
1036
1037 reducedAnswer.add(val.at(indx) * coeff * dV, n);
1038 }
1039
1040 this->giveEdgeDofMapping(mask, iEdge);
1041 answer.assemble(reducedAnswer, mask);
1042 } else {
1043 OOFEM_ERROR("unsupported bc type encountered");
1044 }
1045}
1046
1047
1048void
1049TransportElement :: computeSurfaceBCSubVectorAt(FloatArray &answer, Load *load,
1050 int iSurf, TimeStep *tStep, ValueModeType mode, int indx)
1051{
1052 // OOFEM_ERROR("STOP");
1054 OOFEM_ERROR("no surface load support");
1055 }
1056
1057 BoundaryLoad *surfLoad = dynamic_cast< BoundaryLoad * >(load);
1058 if ( surfLoad ) {
1059 FloatArray reducedAnswer, val, gcoords, n;
1060 IntArray mask;
1061
1062 answer.resize( this->giveNumberOfDofManagers() );
1063 answer.zero();
1064
1065 double coeff = 0.;
1066 int approxOrder = surfLoad->giveApproxOrder();
1067
1068 auto iRule = this->giveBoundarySurfaceIntegrationRule(approxOrder, iSurf);
1069 for ( auto &gp: *iRule ) {
1070 if ( load->giveType() == TransmissionBC ) {
1071 coeff = -1.0;
1072 } else if ( load->giveType() == ConvectionBC ) {
1073 if ( surfLoad->propertyMultExpr.isDefined() ) {
1074 if ( emode != HeatTransferEM && emode != Mass1TransferEM ) {
1076 OOFEM_ERROR("Not implemented for >=2 coupled fields");
1077 }
1078 FloatArray unknowns;
1079 IntArray dofid;
1080 this->computeSurfaceNAt( n, iSurf, gp->giveNaturalCoordinates() );
1081 this->giveElementDofIDMask(dofid);
1082 this->giveSurfaceDofMapping(mask, iSurf);
1083 this->computeBoundaryVectorOf(mask, dofid, VM_TotalIntrinsic, tStep, unknowns);
1084 double value = n.dotProduct(unknowns);
1085 coeff = surfLoad->giveProperty('a', tStep, { { "u", value } });
1086 } else {
1087 coeff = surfLoad->giveProperty('a', tStep);
1088 }
1089 } else if ( load->giveType() == RadiationBC ) {
1090 OOFEM_ERROR("Not implemented, use TransientTransport solver");
1091 //coeff = getRadiativeHeatTranferCoef(surfLoad, tStep);
1092 } else {
1093 OOFEM_ERROR("Unknown load type");
1094 }
1095
1096 this->computeSurfaceNAt( n, iSurf, gp->giveNaturalCoordinates() );
1097 double dV = this->computeSurfaceVolumeAround(gp, iSurf);
1098
1099 FieldPtr tf;
1100 if ( (tf = domain->giveEngngModel()->giveContext()->giveFieldManager()->giveField(FT_TemperatureAmbient)) ) {
1101 //this->computeSurfIpGlobalCoords(gcoords, gp->giveNaturalCoordinates(), iSurf);
1102 this->giveInterpolation()->boundarySurfaceLocal2global(gcoords, iSurf, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this));
1103 tf->evaluateAt(val, gcoords, VM_TotalIntrinsic, tStep);
1104 } else if ( surfLoad->giveFormulationType() == Load :: FT_Entity ) {
1105 surfLoad->computeValueAt(val, tStep, gp->giveNaturalCoordinates(), mode);
1106 } else {
1107 //this->computeSurfIpGlobalCoords(gcoords, gp->giveNaturalCoordinates(), iSurf);
1108 this->giveInterpolation()->boundarySurfaceLocal2global(gcoords, iSurf, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this));
1109
1110 surfLoad->computeValueAt(val, tStep, gcoords, mode);
1111 }
1112
1113 reducedAnswer.add(val.at(indx) * coeff * dV, n);
1114 }
1115
1116 this->giveSurfaceDofMapping(mask, iSurf);
1117 answer.assemble(reducedAnswer, mask);
1118 } else {
1119 OOFEM_ERROR("unsupported bc type encountered");
1120 }
1121}
1122
1123
1124void
1125TransportElement :: computeBCSubMtrxAt(FloatMatrix &answer, TimeStep *tStep, ValueModeType mode, int indx)
1126{
1127 int defined = 0;
1128
1129 answer.resize( this->giveNumberOfDofManagers(), this->giveNumberOfDofManagers() );
1130 answer.zero();
1131
1132 // loop over boundary load array
1133 int nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
1134 for ( int i = 1; i <= nLoads; i++ ) {
1135 int k = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
1136 int id = boundaryLoadArray.at(i * 2);
1137 GeneralBoundaryCondition *load = domain->giveLoad(k);
1138 if ( load->giveType() == ConvectionBC ) {
1139 bcGeomType ltype = load->giveBCGeoType();
1140 if ( ltype == EdgeLoadBGT ) {
1141 BoundaryLoad *edgeLoad = static_cast< BoundaryLoad * >(load);
1142 if ( edgeLoad->isDofExcluded(indx) || !edgeLoad->isImposed(tStep) ) {
1143 continue;
1144 }
1145
1146 defined = 1;
1147 int approxOrder = 2 * this->giveApproxOrder(indx);
1148 int numberOfEdgeIPs = ( int ) ceil( ( approxOrder + 1. ) / 2. );
1149 GaussIntegrationRule iRule(1, this, 1, 1);
1150 iRule.SetUpPointsOnLine(numberOfEdgeIPs, _Unknown);
1151 FloatArray n;
1152 IntArray mask;
1153 FloatMatrix subAnswer;
1154
1155 for ( auto &gp: iRule ) {
1156 this->computeEgdeNAt( n, id, gp->giveNaturalCoordinates() );
1157 double dV = this->computeEdgeVolumeAround(gp, id);
1158 if ( load->giveType() == ConvectionBC ) {
1159 double a;
1160 if( edgeLoad->propertyMultExpr.isDefined() ) {
1161 if ( emode != HeatTransferEM && emode != Mass1TransferEM ) {
1163 OOFEM_ERROR("Not implemented for >=2 coupled fields");
1164 }
1165 FloatArray unknowns;
1166 IntArray dofid;
1167 this->giveElementDofIDMask(dofid);
1168 this->giveEdgeDofMapping(mask, id);
1169 this->computeBoundaryVectorOf(mask, dofid, VM_TotalIntrinsic, tStep, unknowns);
1170 double value = n.dotProduct(unknowns);//unknown in IP
1171 a = edgeLoad->giveProperty('a', tStep, { { "u", value } });
1172 } else {
1173 a = edgeLoad->giveProperty('a', tStep);
1174 }
1175 subAnswer.plusDyadSymmUpper( n, dV * a );
1176 }
1177 }
1178
1179 subAnswer.symmetrized();
1180 this->giveEdgeDofMapping(mask, id);
1181 answer.assemble(subAnswer, mask);
1182 } else if ( ltype == SurfaceLoadBGT ) {
1183 FloatArray n;
1184 IntArray mask;
1185 FloatMatrix subAnswer;
1186
1187 BoundaryLoad *surfLoad = static_cast< BoundaryLoad * >(load);
1188 if ( surfLoad->isDofExcluded(indx) || !surfLoad->isImposed(tStep) ) {
1189 continue;
1190 }
1191
1192 defined = 1;
1193 int approxOrder = 2;
1194 auto iRule = this->giveBoundarySurfaceIntegrationRule(approxOrder, id);
1195 for ( auto &gp: *iRule ) {
1196 this->computeSurfaceNAt( n, id, gp->giveNaturalCoordinates() );
1197 double dV = this->computeSurfaceVolumeAround(gp, id);
1198 if ( load->giveType() == ConvectionBC ) {
1199 double a;
1200 if ( surfLoad->propertyMultExpr.isDefined() ) {
1201 if ( emode != HeatTransferEM && emode != Mass1TransferEM ) {
1203 OOFEM_ERROR("Not implemented for >=2 coupled fields");
1204 }
1205 FloatArray unknowns;
1206 IntArray dofid;
1207 this->giveElementDofIDMask(dofid);
1208 this->giveSurfaceDofMapping(mask, id);
1209 this->computeBoundaryVectorOf(mask, dofid, VM_TotalIntrinsic, tStep, unknowns);
1210 double value = n.dotProduct(unknowns);//unknown in IP
1211 a = surfLoad->giveProperty('a', tStep, { { "u", value } });
1212 } else {
1213 a = surfLoad->giveProperty('a', tStep);
1214 }
1215 subAnswer.plusDyadSymmUpper( n, dV * a );
1216 }
1217 }
1218
1219 subAnswer.symmetrized();
1220 this->giveSurfaceDofMapping(mask, id);
1221 answer.assemble(subAnswer, mask);
1222 } else {
1223 OOFEM_ERROR("unsupported bc type encountered");
1224 }
1225 }
1226 }
1227
1228 // end loop over applied bc
1229
1230 if ( !defined ) {
1231 answer.clear();
1232 }
1233}
1234
1235
1236void
1237TransportElement :: assembleLocalContribution(FloatMatrix &answer, FloatMatrix &src,
1238 int ndofs, int rdof, int cdof)
1239{
1240 int nnodes = this->giveNumberOfDofManagers();
1241
1242 for ( int i = 1; i <= nnodes; i++ ) {
1243 int ti = ( i - 1 ) * ndofs + rdof;
1244 for ( int j = 1; j <= nnodes; j++ ) {
1245 int tj = ( j - 1 ) * ndofs + cdof;
1246 answer.at(ti, tj) += src.at(i, j);
1247 }
1248 }
1249}
1250
1251
1252void
1253TransportElement :: assembleLocalContribution(FloatArray &answer, FloatArray &src,
1254 int ndofs, int rdof)
1255{
1256 int nnodes = this->giveNumberOfDofManagers();
1257
1258 for ( int i = 1; i <= nnodes; i++ ) {
1259 int ti = ( i - 1 ) * ndofs + rdof;
1260 answer.at(ti) += src.at(i);
1261 }
1262}
1263
1264void
1265TransportElement :: computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
1266{
1267 FloatArray u;
1268 FloatMatrix n;
1269 this->computeNmatrixAt(n, lcoords);
1270 this->computeVectorOf(mode, tStep, u);
1271 answer.beProductOf(n, u);
1272}
1273
1274void
1275TransportElement :: computeFlow(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
1276{
1277 FloatArray r, br;
1278 FloatMatrix b, d;
1279 IntArray dofid;
1280
1281 this->giveElementDofIDMask(dofid);
1282 this->computeVectorOf(dofid, VM_TotalIntrinsic, tStep, r);
1284
1285 if ( emode == HeatTransferEM || emode == Mass1TransferEM ) {
1286 this->computeConstitutiveMatrixAt(d, Conductivity_hh, gp, tStep);
1287 br.beProductOf(b, r);
1288 answer.beProductOf(d, br);
1289 } else if ( emode == HeatMass1TransferEM ) {
1290 FloatArray r_h, r_w;
1291 FloatMatrix b_tot, d_hh, d_hw, d_wh, d_ww;
1292 r_h.resize(r.giveSize() / 2);
1293 r_h.zero();
1294 r_w.resize(r.giveSize() / 2);
1295 r_w.zero();
1296
1297 this->computeConstitutiveMatrixAt(d_hh, Conductivity_hh, gp, tStep);
1298 this->computeConstitutiveMatrixAt(d_hw, Conductivity_hw, gp, tStep);
1299 this->computeConstitutiveMatrixAt(d_wh, Conductivity_wh, gp, tStep);
1300 this->computeConstitutiveMatrixAt(d_ww, Conductivity_ww, gp, tStep);
1301
1302 b_tot.resize( 2 * b.giveNumberOfRows(), 2 * b.giveNumberOfColumns() );
1303 b_tot.setSubMatrix(b, 1, 1);
1304 b_tot.setSubMatrix(b, b.giveNumberOfRows() + 1, b.giveNumberOfColumns() + 1);
1305 br.beProductOf(b_tot, r);
1306
1307 d.resize( 2 * d_hh.giveNumberOfRows(), 2 * d_hh.giveNumberOfColumns() );
1308 d.setSubMatrix(d_hh, 1, 1);
1309 d.setSubMatrix(d_hw, 1, d_hh.giveNumberOfColumns() + 1);
1310 d.setSubMatrix(d_wh, d_hh.giveNumberOfRows() + 1, 1);
1311 d.setSubMatrix(d_ww, d_hh.giveNumberOfRows() + 1, d_wh.giveNumberOfColumns() + 1);
1312
1313 answer.beProductOf(d, br);
1314 } else {
1315 OOFEM_ERROR("Unknown element mode");
1316 }
1317
1318 answer.negated();
1319
1320 if ( !isActivated(tStep) ) {
1321 answer.zero();
1322 }
1323}
1324
1325
1326void
1327TransportElement :: updateInternalState(TimeStep *tStep)
1328{
1330 FloatArray stateVector, r;
1331#if 0
1332 FloatArray gradient, flux;
1333 FloatMatrix B;
1334#endif
1335 FloatMatrix n;
1336 IntArray dofid;
1337 TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
1338
1339#ifdef __CEMHYD_MODULE //not very efficient here, looping over all GPs in each call
1340 if (dynamic_cast< CemhydMat * >( mat ) && tStep->isIcApply()){
1341 CemhydMat *cem = dynamic_cast< CemhydMat * >( mat );
1342 cem->initMaterial(this); //create microstructures and statuses on specific GPs
1343 }
1344#endif //__CEMHYD_MODULE
1345
1346 this->giveElementDofIDMask(dofid);
1347 if ( tStep->isIcApply() ) {
1348 this->computeVectorOf(dofid, VM_Total, tStep, r);
1349 } else {
1350 this->computeVectorOf(dofid, VM_TotalIntrinsic, tStep, r);
1351 }
1352 // force updating ip values
1353 for ( auto &iRule: integrationRulesArray ) {
1354 for ( auto &gp: *iRule ) {
1355
1357 this->computeNmatrixAt( n, gp->giveNaturalCoordinates() );
1358 stateVector.beProductOf(n, r);
1359 mat->updateInternalState(stateVector, gp, tStep);
1360
1362#if 0
1363 this->computeGradientMatrixAt( B, gp->giveNaturalCoordinates() );
1364 gradient.beProductOf(B, r);
1365 mat->giveFluxVector(flux, gp, gradient, tStep);
1366#endif
1367 }
1368 }
1369}
1370
1371int
1372TransportElement :: EIPrimaryFieldI_evaluateFieldVectorAt(FloatArray &answer, PrimaryField &pf,
1373 const FloatArray &coords, IntArray &dofId, ValueModeType mode,
1374 TimeStep *tStep)
1375{
1376 int indx;
1377 FloatArray elemvector, lc;
1378 FloatMatrix n;
1379 IntArray elemdofs;
1380 // determine element dof ids
1381 this->giveElementDofIDMask(elemdofs);
1382 // first evaluate element unknown vector
1383 this->computeVectorOf(pf, elemdofs, mode, tStep, elemvector);
1384 // determine corresponding local coordinates
1385 if ( this->computeLocalCoordinates(lc, coords) ) {
1386 // compute interpolation matrix
1387 this->computeNmatrixAt(n, lc);
1388 // compute answer
1389 answer.resize( dofId.giveSize() );
1390 answer.zero();
1391 for ( int i = 1; i <= dofId.giveSize(); i++ ) {
1392 if ( ( indx = elemdofs.findFirstIndexOf( dofId.at(i) ) ) ) {
1393 double sum = 0.0;
1394 for ( int j = 1; j <= elemvector.giveSize(); j++ ) {
1395 sum += n.at(indx, j) * elemvector.at(j);
1396 }
1397
1398 answer.at(i) = sum;
1399 } else {
1400 //_error("EIPrimaryFieldI_evaluateFieldVectorAt: unknown dof id encountered");
1401 answer.at(i) = 0.0;
1402 }
1403 }
1404
1405 return 0; // ok
1406 } else {
1407 OOFEM_ERROR("target point not in receiver volume");
1408 return 1; // failed
1409 }
1410}
1411
1412
1414TransportElement :: giveTransportCrossSection()
1415{
1416 return static_cast< TransportCrossSection* >( this->giveCrossSection() );
1417}
1418
1419
1420Material *
1421TransportElement :: giveMaterial()
1422{
1423 // This method is overloaded in order to have some backwards compatibility
1424 if ( material ) {
1425 return domain->giveMaterial(this->material);
1426 }
1427 return static_cast< TransportCrossSection* >( this->giveCrossSection() )->giveMaterial();
1428}
1429
1430double
1432{
1433 if ( this->vofFunction ) {
1434 return this->giveDomain()->giveFunction(this->vofFunction)->evaluateAtTime(tStep->giveTargetTime());
1435 } else {
1436 return 1.0;
1437 }
1438}
1439
1440#ifdef __OOFEG
1441int
1442TransportElement :: giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode,
1443 int node, TimeStep *tStep)
1444{
1445 Node *n = this->giveNode(node);
1446 if ( type == IST_Temperature ) {
1447 auto dofindx = n->findDofWithDofId(T_f);
1448 if ( dofindx != n->end() ) {
1449 answer.resize(1);
1450 answer.at(1) = (*dofindx)->giveUnknown(VM_Total, tStep);
1451 return 1;
1452 } else {
1453 return 0;
1454 }
1455 } else if ( type == IST_MassConcentration_1 ) {
1456 auto dofindx = n->findDofWithDofId(C_1);
1457 if ( dofindx != n->end() ) {
1458 answer.resize(1);
1459 answer.at(1) = (*dofindx)->giveUnknown(VM_Total, tStep);
1460 return 1;
1461 } else {
1462 return 0;
1463 }
1464 } else {
1465 return Element :: giveInternalStateAtNode(answer, type, mode, node, tStep);
1466 }
1467}
1468
1469#endif
1470} // end namespace oofem
#define N(a, b)
bcType giveType() const override
void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode) override
virtual double giveProperty(int aProperty, TimeStep *tStep, const std ::map< std ::string, FunctionArgument > &valDict) const
virtual double giveTemperOffset(void)
Return temperature offset.
ScalarFunction propertyMultExpr
Expression to multiply all properties.
int giveApproxOrder() override=0
int initMaterial(Element *element) override
Definition cemhydmat.C:352
std::vector< Dof * >::const_iterator findDofWithDofId(DofIDItem dofID) const
Definition dofmanager.C:274
std::vector< Dof * >::iterator end()
Definition dofmanager.h:162
Function * giveFunction(int n)
Definition domain.C:271
virtual void giveElementDofIDMask(IntArray &answer) const
Definition element.h:510
Node * giveNode(int i) const
Definition element.h:629
IntArray boundaryLoadArray
Definition element.h:147
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Definition element.C:1225
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
virtual FEInterpolation * giveInterpolation() const
Definition element.h:648
virtual std::unique_ptr< IntegrationRule > giveBoundaryEdgeIntegrationRule(int order, int boundary)
Definition element.C:911
virtual int computeNumberOfDofs()
Definition element.h:436
void initializeFrom(InputRecord &ir, int priority) override
Definition element.C:687
virtual std::unique_ptr< IntegrationRule > giveBoundarySurfaceIntegrationRule(int order, int boundary)
Definition element.C:917
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
IntArray * giveBoundaryLoadArray()
Returns array containing load numbers of boundary loads acting on element.
Definition element.C:420
void postInitialize() override
Performs post initialization steps.
Definition element.C:797
virtual IntegrationRule * giveIntegrationRule(int i)
Definition element.h:899
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int material
Number of associated material.
Definition element.h:140
virtual int giveNumberOfDofManagers() const
Definition element.h:695
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
Element(int n, Domain *aDomain)
Definition element.C:86
virtual int testElementExtension(ElementExtension ext)
Definition element.h:910
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Definition element.C:1240
virtual double computeVolumeAround(GaussPoint *gp)
Definition element.h:530
void computeBoundaryVectorOf(const IntArray &bNodes, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false)
Definition element.C:163
IntArray bodyLoadArray
Definition element.h:147
virtual Element_Geometry_Type giveGeometryType() const =0
virtual IntArray boundaryEdgeGiveNodes(int boundary, const Element_Geometry_Type, bool includeHierarchical=false) const =0
virtual IntArray boundaryGiveNodes(int boundary, const Element_Geometry_Type) const =0
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void boundaryEdgeEvalN(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual std::unique_ptr< IntegrationRule > giveIntegrationRule(int order, const Element_Geometry_Type) const
Definition feinterpol.C:90
int giveInterpolationOrder() const
Definition feinterpol.h:199
virtual void boundaryLocal2Global(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
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
int number
Component number.
Definition femcmpnn.h:77
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
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void power(const double exponent)
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
bool isEmpty() const
Returns true if receiver is empty.
Definition floatarray.h:265
void add(const FloatArray &src)
Definition floatarray.C:218
void subtract(const FloatArray &src)
Definition floatarray.C:320
void times(double s)
Definition floatarray.C:834
bool isNotEmpty() const
Returns true if receiver is not empty.
Definition floatarray.h:263
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
void add(const FloatMatrix &a)
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()
int giveNumberOfColumns() const
Returns number of columns of receiver.
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void beNMatrixOf(const FloatArray &n, int nsd)
void setSubMatrix(const FloatMatrix &src, int sr, int sc)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
bool isNotEmpty() const
Tests for empty matrix.
void plusDyadSymmUpper(const FloatArray &a, double dV)
int giveNumberOfRows() const
Returns number of rows of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
double at(std::size_t i, std::size_t j) const
virtual double evaluateAtTime(double t)
Definition function.C:78
int SetUpPointsOnLine(int nPoints, MaterialMode mode) override
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
virtual bcGeomType giveBCGeoType() const
virtual bool isImposed(TimeStep *tStep)
int findFirstIndexOf(int value) const
Definition intarray.C:280
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)=0
virtual int giveApproxOrder()
Definition load.h:174
int isDofExcluded(int index)
Definition load.C:136
virtual FormulationType giveFormulationType()
Definition load.h:170
double giveTargetTime()
Returns target time.
Definition timestep.h:164
bool isIcApply()
Definition timestep.C:154
static ParamKey IPK_TransportElement_vof_function
virtual void computeEgdeNAt(FloatArray &answer, int iEdge, const FloatArray &lcoord)
virtual void computeInertiaForcesVector(FloatArray &answer, TimeStep *tStep)
void computeSurfaceBCSubVectorAt(FloatArray &answer, Load *load, int iSurf, TimeStep *tStep, ValueModeType mode, int indx)
void computeBCSubVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode, int indx)
static const double stefanBoltzmann
Stefan–Boltzmann constant W/m2/K4.
virtual void giveSurfaceDofMapping(IntArray &mask, int iSurf)
void computeBCSubMtrxAt(FloatMatrix &answer, TimeStep *tStep, ValueModeType mode, int indx)
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
virtual void computeConductivityMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
void computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge, TimeStep *tStep, ValueModeType mode, int indx)
virtual void computeSurfaceNAt(FloatArray &answer, int iSurf, const FloatArray &lcoord)
virtual double computeVof(TimeStep *tStep)
Material * giveMaterial() override
virtual void computeBmatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
virtual double giveThicknessAt(const FloatArray &gcoords)
int vofFunction
Fuction determining the relative volume of reference material in element.
virtual void computeNmatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
virtual int giveApproxOrder(int unknownIndx)
virtual void computeIntSourceLHSSubMatrix(FloatMatrix &answer, MatResponseMode rmode, int iri, TimeStep *tStep)
void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep) override
virtual void computeConductivitySubMatrix(FloatMatrix &answer, int iri, MatResponseMode rmode, TimeStep *tStep)
void computeBodyBCSubVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, int indx)
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
virtual void computeNAt(FloatArray &answer, const FloatArray &lcoord)
virtual void giveEdgeDofMapping(IntArray &mask, int iEdge)
virtual void computeCapacityMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual void computeCapacitySubMatrix(FloatMatrix &answer, MatResponseMode rmode, int iri, TimeStep *tStep)
virtual void computeInternalSourceRhsSubVectorAt(FloatArray &answer, TimeStep *, ValueModeType mode, int indx)
virtual void computeGradientMatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
virtual void computeLumpedCapacityVector(FloatArray &answer, TimeStep *tStep)
virtual void computeBCMtrxAt(FloatMatrix &answer, TimeStep *tStep, ValueModeType mode)
void assembleLocalContribution(FloatMatrix &answer, FloatMatrix &src, int ndofs, int rdof, int cdof)
virtual void computeBCVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
virtual void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
virtual void updateInternalState(const FloatArray &state, GaussPoint *gp, TimeStep *tStep)
virtual void computeInternalSourceVector(FloatArray &val, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const
virtual double giveCharacteristicValue(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override=0
virtual void giveFluxVector(FloatArray &answer, GaussPoint *gp, const FloatArray &grad, const FloatArray &field, TimeStep *tStep) const
virtual bool hasInternalSource() const
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
bcGeomType
Type representing the geometric character of loading.
Definition bcgeomtype.h:40
@ SurfaceLoadBGT
Distributed surface load.
Definition bcgeomtype.h:45
@ EdgeLoadBGT
Distributed edge load.
Definition bcgeomtype.h:44
const char * __CharTypeToString(CharType _value)
Definition cltypes.C:338
InternalStateMode
Determines the mode of internal variable.
double sum(const FloatArray &x)
@ RadiationBC
Stefan-Boltzmann law.
Definition bctype.h:48
@ TransmissionBC
Neumann type (prescribed flux).
Definition bctype.h:43
@ ConvectionBC
Newton type - transfer coefficient.
Definition bctype.h:44
@ Element_SurfaceLoadSupport
Element extension for surface loads.
std::shared_ptr< Field > FieldPtr
Definition field.h:78
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)

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