OOFEM 3.0
Loading...
Searching...
No Matches
macrolspace.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
39#include "fei3dhexalin.h"
40#include "constantfunction.h"
41#include "domain.h"
42#include "classfactory.h"
43#include "dynamicinputrecord.h"
44#include "node.h"
45#include "parametermanager.h"
46#include "paramkey.h"
47
48namespace oofem {
53
54//derived from linear brick element
55MacroLSpace :: MacroLSpace(int n, Domain *aDomain) : LSpace(n, aDomain)
56{
57 this->firstCall = true;
58 microMaterial = NULL;
59 microDomain = NULL;
60 microEngngModel = NULL;
61 this->iteration = 1;
62 this->lastStiffMatrixTimeStep = NULL;
63}
64
65MacroLSpace :: ~MacroLSpace() { }
66
67
68void MacroLSpace :: initializeFrom(InputRecord &ir, int priority)
69{
70 LSpace :: initializeFrom(ir, priority);
71 ParameterManager &ppm = domain->elementPPM;
74
75
76#if 0 // obsolete
79 if ( fopen(this->stiffMatrxFileName, "r") != NULL ) { //if the file exist
80 stiffMatrxFile = fopen(this->stiffMatrxFileName, "r");
82 } else { //or create a new one
83 if ( ( stiffMatrxFile = fopen(this->stiffMatrxFileName, "w") ) == NULL ) {
84 OOFEM_ERROR("Can not create a new file %s\n", this->stiffMatrxFileName);
85 }
87 }
88 }
89#endif
90}
91
92void MacroLSpace :: postInitialize()
93{
94 LSpace :: postInitialize();
95 ParameterManager &ppm = domain->elementPPM;
98 if ( this->microMasterNodes.giveSize() != 8 ) {
99 throw ComponentInputException(IPK_MacroLSpace_microMasterNodes.getName(), ComponentInputException::ComponentType::ctElement, this->number, "Need 8 master nodes from the microproblem defined on macroLspace element");
100 }
101 microBoundaryDofManager.resize( 3 * microBoundaryNodes.giveSize() );
102}
103
104
105/*Stiffness matrix is taken from the microproblem. No GPs are presented here on macroscale.
106 * Stiffness matrix (24,24) goes in order node 1 (u,v,w) to node 2 (u,v,w) ... 8 (u,v,w). Displacements are in global coordinates.
107 */
108void MacroLSpace :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
109{
110 //MatResponseMode rMode specifies tangent, secant, or initial matrix
111 if ( !this->isActivated(tStep) ) {
112 return;
113 }
114
115
116 //called the first time and initiates the microproblem
117 if ( this->firstCall ) {
118 this->microMaterial = static_cast< MicroMaterial * >( this->giveMaterial() );
119 if ( this->microMaterial->microMatIsUsed == true ) {
120 OOFEM_ERROR("Micromaterial is already used on another element. Only one micromaterial can be assigned to one macro element");
121 }
122
123 this->microDomain = this->microMaterial->problemMicro->giveDomain(1); //from engngm.h
124 this->microEngngModel = this->microDomain->giveEngngModel();
125 this->microEngngModel->setProblemScale(microScale); //set microScale attribute
126 this->microEngngModel->checkProblemConsistency();
127 this->microMaterial->init(); //from UnknownNumberingScheme()
128 this->microMaterial->setMacroProperties(this->giveDomain(), this, this->microMasterNodes, this->microBoundaryNodes);
129 this->firstCall = false;
130 }
131
132 //call microproblem
133 //activeMStep = microMat->problemMicro->giveMetaStep(1);//->setNumberOfSteps(1);
134 //activeMStep->giveMetaStepNumber();
135 //microEngngModel->timer.startTimer(EngngModelTimer :: EMTT_AnalysisTimer);
136 //microproblem must have the same actual time and zero time increment
137
138 this->microEngngModel->giveCurrentStep()->setTargetTime( tStep->giveTargetTime() ); //adjust total time
139 this->microEngngModel->giveCurrentStep()->setIntrinsicTime( tStep->giveIntrinsicTime() ); //adjust intrinsic time
140 this->microEngngModel->giveCurrentStep()->setTimeIncrement(0.); //no time increment
141 this->microEngngModel->initMetaStepAttributes( microEngngModel->giveCurrentMetaStep() ); //updates numerical method
142
143 OOFEM_LOG_INFO( "\n** Assembling %s stiffness matrix of microproblem on macroElement %d, micTimeStep %d, micTime %f\n", __MatResponseModeToString(rMode), this->giveNumber(), this->microEngngModel->giveCurrentStep()->giveNumber(), this->microEngngModel->giveCurrentStep()->giveTargetTime() );
144
145 //this->microEngngModel->solveYourselfAt( microEngngModel->giveCurrentStep() );
146 //this->microEngngModel->terminate( microEngngModel->giveCurrentStep() );
147
148 if ( tStep != this->lastStiffMatrixTimeStep ) {
149 this->microMaterial->giveMacroStiffnessMatrix(answer, this->microEngngModel->giveCurrentStep(), rMode, this->microMasterNodes, this->microBoundaryNodes);
150 this->stiffMatrix = answer;
151 this->lastStiffMatrixTimeStep = tStep;
152 OOFEM_LOG_INFO("** Assembled now\n\n");
153 } else {
154 answer = this->stiffMatrix;
155 OOFEM_LOG_INFO("** Assembled previously in this time step\n\n");
156 }
157}
158
159
160//assign values to DOF on the boundary according to definition on macrolspace and actual displacement stage
161void MacroLSpace :: changeMicroBoundaryConditions(TimeStep *tStep)
162{
163 DynamicInputRecord ir_func, ir_bc;
164 FloatArray n(8), answer, localCoords;
165 double displ;
166 FloatArray displ_x(8), displ_y(8), displ_z(8);
167 int counter;
168
169 //dofManArray has the node order as specified in input file
170 for ( int i = 1; i <= this->giveNumberOfNodes(); i++ ) { //8 nodes
171 //global displacements
172 displ_x.at(i) = this->giveNode(i)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
173 displ_y.at(i) = this->giveNode(i)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
174 displ_z.at(i) = this->giveNode(i)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
175 }
176
177 //overrides the first load-time function to be a constant
178 if ( !microDomain->giveNumberOfFunctions() ) {
179 microDomain->resizeFunctions(1);
180 }
181
182 ir_func.setRecordKeywordField("constantfunction", 1);
183 ir_func.setField(1.0, _IFT_ConstantFunction_f);
184 auto timeFunct = classFactory.createFunction("constantfunction", 1, microDomain);
185 if ( !timeFunct ) {
186 OOFEM_ERROR("Couldn't create constant time function");
187 }
188 timeFunct->initializeFrom(ir_func);
189 microDomain->setFunction(1, std::move(timeFunct));
190
191
192 /*assign to each boundary node the form "bc 3 # # #", set 0s on free nodes
193 * modify bcList = corresponds to "nbc"
194 */
195 microDomain->clearBoundaryConditions(); //from domain.C
196 microDomain->resizeBoundaryConditions( 3 * microBoundaryNodes.giveSize() ); //from domain.C
197
198 counter = 1;
199 for ( auto &DofMan : microDomain->giveDofManagers() ) { //go through all nodes on microDomain
200 if ( microBoundaryNodes.contains( DofMan->giveGlobalNumber() ) ) { //if the node number is on boundary
201 this->evalInterpolation( n, microMaterial->microMasterCoords, DofMan->giveCoordinates() );
202 //n.printYourself();
203
204 for ( Dof *dof: *DofMan ) {
205 this->microBoundaryDofManager.at(counter) = DofMan->giveGlobalNumber();
206 dof->setBcId(counter);
207 DofIDItem id = dof->giveDofID();
208 displ = n.dotProduct( id == D_u ? displ_x : ( id == D_v ? displ_y : displ_z ) );
209 ir_bc.setRecordKeywordField("boundarycondition", counter);
212 auto bc = classFactory.createBoundaryCondition("boundarycondition", counter, microDomain);
213 if ( !bc ) {
214 OOFEM_ERROR("Couldn't create boundary condition.");
215 }
216 bc->initializeFrom(ir_bc);
217 microDomain->setBoundaryCondition(counter, std::move(bc));
218 counter++;
219 }
220 } else {
221 for ( auto &dof: *DofMan ) {
222 dof->setBcId(0);
223 }
224 }
225 }
226}
227
228//obtain nodal forces from underlying microScale
229//node numbering on element is in the same order as in the input file
230//useUpdatedGpRecord=1 is used for printing of reactions
231void MacroLSpace :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
232{
233 //StructuralEngngModel *microStructuralEngngModel;
234 FloatArray reactions, localCoords, n(8);
235 DofManager *DofMan;
236 double reactionForce;
237
238 this->microEngngModel->giveCurrentStep()->setTargetTime( tStep->giveTargetTime() ); //adjust total time
239 this->microEngngModel->giveCurrentStep()->setIntrinsicTime( tStep->giveIntrinsicTime() ); //adjust total time
240 this->microEngngModel->giveCurrentStep()->setTimeIncrement(0.); //no time increment
241
242 //OOFEM_LOG_INFO("*** useUpdatedGpRecord %d\n", useUpdatedGpRecord);
243
244 if ( useUpdatedGpRecord ) { //printing of data
245 answer = this->internalMacroForcesVector;
246 } else {
247 OOFEM_LOG_INFO( "\n*** Solving reactions of macroElement %d, micTimeStep %d, macIteration %d, micTime %f\n", this->giveNumber(), this->microEngngModel->giveCurrentStep()->giveNumber(), this->iteration, this->microEngngModel->giveCurrentStep()->giveTargetTime() );
248
249 this->iteration++;
251
252 this->microEngngModel->solveYourselfAt( this->microEngngModel->giveCurrentStep() );
253 this->microEngngModel->updateYourself( this->microEngngModel->giveCurrentStep() );
254 //this->microEngngModel->terminate( this->microEngngModel->giveCurrentStep() );
255 StructuralEngngModel *microStructuralEngngModel = dynamic_cast< StructuralEngngModel * >(this->microEngngModel);
256
257 //reaction vector contains contributions from unknownNumberingScheme
258 microStructuralEngngModel->computeReaction(reactions, this->microEngngModel->giveCurrentStep(), 1);
259 //reactions.printYourself();
260 answer.resize(24);
261 answer.zero();
262 /*for ( i = 1; i <= this->giveNumberOfNodes(); i++ ) { //8 nodes
263 * DofMan = microDomain->giveDofManager(i);
264 */
265
266 for ( int i = 1; i <= this->microBoundaryDofManager.giveSize() / 3; i++ ) { //Number of DoFManagers stored in triplets
267 DofMan = microDomain->giveDofManager( this->microBoundaryDofManager.at(3 * i - 2) );
268 this->evalInterpolation( n, microMaterial->microMasterCoords, DofMan->giveCoordinates() );
269 for ( int j = 1; j <= DofMan->giveNumberOfDofs(); j++ ) { //3
270 reactionForce = reactions.at(3 * i + j - 3);
271 for ( int k = 1; k <= 8; k++ ) {
272 answer.at(3 * k + j - 3) += reactionForce * n.at(k);
273 }
274 }
275 }
276
277 this->internalMacroForcesVector = answer;
278 OOFEM_LOG_INFO("*** Reactions done\n");
279 }
280
281 //answer.printYourself();
282 //OOFEM_ERROR("STOP");
283}
284
285void MacroLSpace :: evalInterpolation(FloatArray &answer, const std::vector< FloatArray > &coords, const FloatArray &gcoords)
286{
287 FloatArray localCoords;
288
289 //this->interpolation.global2local(localCoords, coords, gcoords, 0.0);//returns even outside the element boundaries
290 //this->interpolation.evalN(answer, localCoords, 0.0);
291 this->interpolation.global2local( localCoords, gcoords, FEIVertexListGeometryWrapper(coords, this->interpolation.giveGeometryType()) ); //returns even outside the element boundaries
292 this->interpolation.evalN( answer, localCoords, FEIVertexListGeometryWrapper(coords, this->interpolation.giveGeometryType()) );
293}
294
295
296// Updates the receiver at end of time step.
297void MacroLSpace :: updateYourself(TimeStep *tStep)
298{
299 FloatArray answer;
300
301 for ( auto &iRule: integrationRulesArray ) {
302 iRule->updateYourself(tStep);
303 }
304
305 OOFEM_LOG_INFO("*** Updating macroelement\n");
306 //set number of timestep so the microproblem counts from 1 to nsteps in each iteration but not totally
307 //this->microEngngModel->giveCurrentStep()->setNumber(1);
308 //recalculate microproblem
309 //this->giveInternalForcesVector(answer, tStep, 0);
310
311 this->microEngngModel->terminate( this->microEngngModel->giveCurrentStep() ); //perform output, VTK
312 this->iteration = 1;
313 this->microEngngModel->giveNextStep(); //time step used in ouput file name
314}
315} // end namespace oofem
#define _IFT_BoundaryCondition_PrescribedValue
[rn,optional] Prescribed value of all DOFs
#define REGISTER_Element(class)
int giveNumberOfDofs() const
Definition dofmanager.C:287
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
void setRecordKeywordField(std ::string keyword, int number)
void setField(int item, InputFieldType id)
Node * giveNode(int i) const
Definition element.h:629
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
virtual int giveNumberOfNodes() const
Definition element.h:703
virtual Material * giveMaterial()
Definition element.C:523
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
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
int giveNumber() const
Definition femcmpnn.h:104
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
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
LSpace(int n, Domain *d)
Definition lspace.C:62
static FEI3dHexaLin interpolation
Definition lspace.h:67
virtual void evalInterpolation(FloatArray &answer, const std::vector< FloatArray > &coords, const FloatArray &gcoords)
IntArray microBoundaryNodes
Definition macrolspace.h:99
IntArray microMasterNodes
Array containing the node mapping from microscale (which microMasterNodes corresponds to which macroN...
Definition macrolspace.h:98
IntArray microBoundaryDofManager
Stores node number on the boundary in the triplets.
MicroMaterial * microMaterial
FloatMatrix stiffMatrix
EngngModel * microEngngModel
FloatArray internalMacroForcesVector
Array containg the force vector from nodes (if condensation is skipped, use this vector).
static ParamKey IPK_MacroLSpace_stiffMatrixFileName
Definition macrolspace.h:65
TimeStep * lastStiffMatrixTimeStep
Last time step when stiffness matrix was assembled.
static ParamKey IPK_MacroLSpace_microMasterNodes
Definition macrolspace.h:63
virtual void changeMicroBoundaryConditions(TimeStep *tStep)
Related to setting the boundary conditions of micro problem.
int stiffMatrxFileNoneReadingWriting
Process with external file for the storage of stiffness matrix 0-None, 1-read, 2-write.
static ParamKey IPK_MacroLSpace_microBoundaryNodes
Definition macrolspace.h:64
int iteration
Information of iteration number.
void computeReaction(FloatArray &answer, TimeStep *tStep, int di)
double giveTargetTime()
Returns target time.
Definition timestep.h:164
double giveIntrinsicTime()
Returns intrinsic time, e.g. time in which constitutive model is evaluated.
Definition timestep.h:166
#define _IFT_ConstantFunction_f
#define OOFEM_ERROR(...)
Definition error.h:79
#define _IFT_GeneralBoundaryCondition_timeFunct
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
#define _IFT_MacroLspace_stiffMatrxFileName
Definition macrolspace.h:45
ClassFactory & classFactory
@ microScale
Definition problemmode.h:47
#define PM_ELEMENT_ERROR_IFNOTSET(_pm, _componentnum, _paramkey)
#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