OOFEM 3.0
Loading...
Searching...
No Matches
tr1darcy.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
37#include "fei2dtrlin.h"
38#include "floatarrayf.h"
39#include "floatmatrixf.h"
41#include "gausspoint.h"
42#include "bcgeomtype.h"
45#include "load.h"
46#include "boundaryload.h"
47#include "mathfem.h"
48#include "crosssection.h"
49#include "classfactory.h"
50
51#include "iostream"
52
53namespace oofem {
55
56FEI2dTrLin Tr1Darcy :: interpolation_lin(1, 2);
57
58Tr1Darcy :: Tr1Darcy(int n, Domain *aDomain) : TransportElement(n, aDomain)
59{
61 this->numberOfGaussPoints = 1;
62}
63
64void Tr1Darcy :: initializeFrom(InputRecord &ir, int priority)
65{
66 TransportElement :: initializeFrom(ir, priority);
67}
68
70Tr1Darcy :: giveInterpolation() const
71{
72 return & interpolation_lin;
73}
74
75void Tr1Darcy :: computeGaussPoints()
76{
77 if ( integrationRulesArray.size() == 0 ) {
78 integrationRulesArray.resize( 1 );
79 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 3);
80 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
81 }
82}
83
84void Tr1Darcy :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
85{
86 auto mat = static_cast< TransportMaterial * >( this->giveMaterial() );
87
89 for ( auto &gp: *integrationRulesArray [ 0 ] ) {
90 // auto [detJ, B] = evaldNdx(lcoords, FEIElementGeometryWrapper(this));
91 auto x = this->interpolation_lin.evaldNdx(FEIElementGeometryWrapper(this));
92 auto detJ = x.first;
93 auto B = x.second;
94
95 auto D = mat->computeTangent2D(mode, gp, tStep);
96 Ke += detJ * gp->giveWeight() * Tdot(B, dot(D, B));
97 }
98 answer = Ke;
99}
100
101void Tr1Darcy :: giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode, TimeStep *tStep)
102{
103 if ( mtrx == ExternalForcesVector ) {
104 this->computeExternalForcesVector(answer, tStep, mode);
105 } else if ( mtrx == InternalForcesVector ) {
106 this->computeInternalForcesVector(answer, tStep);
107 } else {
108 OOFEM_ERROR("Unknown Type of characteristic mtrx.");
109 }
110}
111
112void Tr1Darcy :: computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
113{
114 auto mat = static_cast< TransportMaterial * >( this->giveMaterial() );
115
116 FloatArray a_tmp;
117 this->computeVectorOf(VM_Total, tStep, a_tmp);
118 FloatArrayF<3> a = a_tmp;
119
121 for ( auto &gp: *integrationRulesArray [ 0 ] ) {
122 const FloatArrayF<2> lcoords = gp->giveNaturalCoordinates();
123
124 // auto [detJ, B] = evaldNdx(lcoords, FEIElementGeometryWrapper(this));
125 auto x = this->interpolation_lin.evaldNdx(FEIElementGeometryWrapper(this));
126 auto detJ = x.first;
127 auto B = x.second;
128 auto n = this->interpolation_lin.evalN(lcoords);
129
130 auto p = dot(n, a);
131 auto gradp = dot(B, a);
132 auto w = mat->computeFlux2D(gradp, p, gp, tStep);
133
134 fe += -gp->giveWeight() * detJ * Tdot(B, w);
135 }
136 answer = fe;
137}
138
139
140void Tr1Darcy :: computeExternalForcesVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
141{
142 // TODO: Implement support for body forces
143
144 FloatArray vec;
145
146 answer.resize(3);
147 answer.zero();
148
149 // Compute characteristic vector for Neumann boundary conditions.
150 int nLoads = boundaryLoadArray.giveSize() / 2;
151 for ( int i = 1; i <= nLoads; i++ ) { // For each Neumann boundary condition ....
152 int load_number = boundaryLoadArray.at(2 * i - 1);
153 int load_id = boundaryLoadArray.at(2 * i);
154 auto load = domain->giveLoad(load_number);
155 auto ltype = load->giveBCGeoType();
156
157 if ( ltype == EdgeLoadBGT ) {
158 this->computeEdgeBCSubVectorAt(vec, load, load_id, tStep, mode, 0);
159 }
160
161 answer.add(vec);
162 }
163
164 answer.negated();
165}
166
167void Tr1Darcy :: computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge, TimeStep *tStep, ValueModeType mode, int indx)
168{
169 /*
170 * Given the load *load, return it's contribution.
171 *
172 */
173
174 answer.resize(3);
175 answer.zero();
176
177 if ( load->giveType() == TransmissionBC ) { // Neumann boundary conditions (traction)
178 auto boundaryLoad = static_cast< BoundaryLoad * >(load);
179
180 int numberOfEdgeIPs;
181 numberOfEdgeIPs = ( int ) ceil( ( boundaryLoad->giveApproxOrder() + 1. ) / 2. ) * 2;
182
183 GaussIntegrationRule iRule(1, this, 1, 1);
184 FloatArray N, loadValue, reducedAnswer(3);
185
186 iRule.SetUpPointsOnLine(numberOfEdgeIPs, _Unknown);
187
188 for ( auto &gp: iRule ) {
189 const FloatArray &lcoords = gp->giveNaturalCoordinates();
190 this->interpolation_lin.edgeEvalN( N, iEdge, lcoords, FEIElementGeometryWrapper(this) );
191 double dV = this->computeEdgeVolumeAround(gp, iEdge);
192
193 if ( boundaryLoad->giveFormulationType() == Load :: FT_Entity ) { // Edge load in xi-eta system
194 boundaryLoad->computeValueAt(loadValue, tStep, lcoords, mode);
195 } else { // Edge load in x-y system
196 FloatArray gcoords;
197 this->interpolation_lin.edgeLocal2global( gcoords, iEdge, lcoords, FEIElementGeometryWrapper(this) );
198 boundaryLoad->computeValueAt(loadValue, tStep, gcoords, mode);
199 }
200
201 reducedAnswer.add(loadValue.at(1) * dV, N);
202 }
203
204 const auto &mask = this->interpolation_lin.computeLocalEdgeMapping(iEdge);
205 answer.assemble(reducedAnswer, mask);
206 }
207}
208
209double Tr1Darcy :: giveThicknessAt(const FloatArray &gcoords)
210{
211 return this->giveCrossSection()->give(CS_Thickness, gcoords, this, false);
212}
213
214
215double Tr1Darcy :: computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
216{
217 double thickness = 1;
218 double detJ = fabs( this->interpolation_lin.edgeGiveTransformationJacobian( iEdge, gp->giveSubPatchCoordinates(), FEIElementGeometryWrapper(this) ) );
219 return detJ *thickness *gp->giveWeight();
220}
221
222void Tr1Darcy :: giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep)
223{
224 if ( mtrx == ConductivityMatrix || mtrx == TangentStiffnessMatrix ) {
225 this->computeStiffnessMatrix(answer, TangentStiffness, tStep);
226 } else {
227 OOFEM_ERROR("Unknown Type of characteristic mtrx %d", mtrx);
228 }
229}
230
231void Tr1Darcy :: giveDofManDofIDMask(int inode, IntArray &answer) const
232{
233 answer = {P_f};
234}
235
236void
237Tr1Darcy :: NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node,
238 InternalStateType type, TimeStep *tStep)
239{
241 auto gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
242 this->giveCrossSection()->giveIPValue(answer, gp, type, tStep);
243}
244
245Interface *
246Tr1Darcy :: giveInterface(InterfaceType interface)
247{
248 if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
249 return static_cast< NodalAveragingRecoveryModelInterface * >(this);
250 }
251
252 return nullptr;
253}
254
255int
256Tr1Darcy :: computeNumberOfDofs()
257{
258 return 3;
259}
260}
#define N(a, b)
#define REGISTER_Element(class)
IntArray boundaryLoadArray
Definition element.h:147
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
CrossSection * giveCrossSection()
Definition element.C:534
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
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 zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void add(const FloatArray &src)
Definition floatarray.C:218
int SetUpPointsOnLine(int nPoints, MaterialMode mode) override
const FloatArray & giveSubPatchCoordinates() const
Returns local sub-patch coordinates of the receiver.
Definition gausspoint.h:142
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep) override
Definition tr1darcy.C:112
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Definition tr1darcy.C:84
void computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge, TimeStep *tStep, ValueModeType mode, int indx)
Definition tr1darcy.C:167
double computeEdgeVolumeAround(GaussPoint *gp, int iEdge) override
Definition tr1darcy.C:215
void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode) override
Definition tr1darcy.C:140
static FEI2dTrLin interpolation_lin
Definition tr1darcy.h:54
TransportElement(int n, Domain *d, ElementMode em=HeatTransferEM)
Material * giveMaterial() override
#define OOFEM_ERROR(...)
Definition error.h:79
@ CS_Thickness
Thickness.
FloatMatrixF< N, P > Tdot(const FloatMatrixF< M, N > &a, const FloatMatrixF< M, P > &b)
Computes .
@ EdgeLoadBGT
Distributed edge load.
Definition bcgeomtype.h:44
double dot(const FloatArray &x, const FloatArray &y)
@ NodalAveragingRecoveryModelInterfaceType
@ TransmissionBC
Neumann type (prescribed flux).
Definition bctype.h:43

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