OOFEM 3.0
Loading...
Searching...
No Matches
fluidstructureproblem.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 "engngm.h"
37#include "timestep.h"
38#include "function.h"
39#include "metastep.h"
40#include "exportmodulemanager.h"
41#include "mathfem.h"
42#include "oofemtxtdatareader.h"
43#include "util.h"
44#include "verbose.h"
45#include "classfactory.h"
46
47#include "pfem.h"
51
52#include <stdlib.h>
53
54#ifdef __OOFEG
55 #include "oofeggraphiccontext.h"
56#endif
57
58namespace oofem {
60
61FluidStructureProblem :: FluidStructureProblem(int i, EngngModel *_master) : StaggeredProblem(i, _master)
62{
63 ndomains = 1; // domain is needed to store the time step ltf
64
66 this->setRenumberFlag();
67}
68
69FluidStructureProblem :: ~FluidStructureProblem()
70{}
71
72
73void
74FluidStructureProblem :: initializeFrom(InputRecord &ir)
75{
76 StaggeredProblem :: initializeFrom(ir);
77
78 maxiter = 50;
80
81 rtolv = 1.e-3;
83
84 rtolp = 1.e-3;
86}
87
88void
89FluidStructureProblem :: initializeYourself(TimeStep *tStep)
90{
91 for ( auto &emodel: emodelList ) {
92 emodel->initializeYourself(tStep);
93
94 DIIDynamic *dynamicProblem = dynamic_cast< DIIDynamic * >( emodel.get() );
95 if ( dynamicProblem ) {
96 this->giveCurrentStep()->setTimeDiscretization( dynamicProblem->giveInitialTimeDiscretization() );
97 }
98 NonLinearDynamic *nonlinearDynamicProblem = dynamic_cast< NonLinearDynamic * >( emodel.get() );
99 if ( nonlinearDynamicProblem ) {
100 this->giveCurrentStep()->setTimeDiscretization( nonlinearDynamicProblem->giveInitialTimeDiscretization() );
101 }
102 }
103}
104
105void
106FluidStructureProblem :: solveYourselfAt(TimeStep *stepN)
107{
108 PFEM *pfemProblem = dynamic_cast< PFEM * >( this->giveSlaveProblem(1) );
109 Domain *pfemDomain = pfemProblem->giveDomain(1);
110
111#ifdef VERBOSE
112 OOFEM_LOG_RELEVANT( "Solving [step number %5d, time %e]\n", stepN->giveNumber(), stepN->giveTargetTime() );
113#endif
114
115 FloatArray currentInteractionParticlesVelocities;
116 FloatArray currentInteractionParticlesPressures;
117 FloatArray previousInteractionParticlesVelocities;
118 FloatArray previousInteractionParticlesPressures;
119
120 if ( stepN->isTheFirstStep() ) {
121 int ndman = pfemDomain->giveNumberOfDofManagers();
122 for ( int i = 1; i <= ndman; i++ ) {
123 if ( dynamic_cast< InteractionPFEMParticle * >( pfemDomain->giveDofManager(i) ) ) {
124 interactionParticles.followedBy(i, 10);
125 }
126 }
127 }
128
129 currentInteractionParticlesVelocities.resize( 2 * interactionParticles.giveSize() );
130 previousInteractionParticlesVelocities.resize( 2 * interactionParticles.giveSize() );
131 currentInteractionParticlesPressures.resize( interactionParticles.giveSize() );
132 previousInteractionParticlesPressures.resize( interactionParticles.giveSize() );
133
134 double velocityDifference = 1.0;
135 double pressureDifference = 1.0;
136
137 iterationNumber = 0;
138 while ( ( ( velocityDifference > rtolv ) || ( pressureDifference > rtolp ) ) && iterationNumber < maxiter ) {
139 previousInteractionParticlesPressures = currentInteractionParticlesPressures;
140 previousInteractionParticlesVelocities = currentInteractionParticlesVelocities;
141
142 for ( auto &emodel: emodelList ) {
143 emodel->solveYourselfAt(stepN);
144 }
145
146 for ( int i = 1; i <= interactionParticles.giveSize(); i++ ) {
147 currentInteractionParticlesPressures.at(i) = pfemProblem->giveUnknownComponent( VM_Total, stepN, pfemDomain, pfemDomain->giveDofManager( interactionParticles.at(i) )->giveDofWithID(P_f) );
148 InteractionPFEMParticle *interactionParticle = dynamic_cast< InteractionPFEMParticle * >( pfemDomain->giveDofManager( interactionParticles.at(i) ) );
149 FloatArray velocities;
150 interactionParticle->giveCoupledVelocities(velocities, stepN);
151 currentInteractionParticlesVelocities.at(2 * ( i - 1 ) + 1) = velocities.at(1);
152 currentInteractionParticlesVelocities.at(2 * ( i - 1 ) + 2) = velocities.at(2);
153 }
154
155 pressureDifference = 0.0;
156 velocityDifference = 0.0;
157 for ( int i = 1; i <= currentInteractionParticlesPressures.giveSize(); i++ ) {
158 pressureDifference += ( currentInteractionParticlesPressures.at(i) - previousInteractionParticlesPressures.at(i) ) * ( currentInteractionParticlesPressures.at(i) - previousInteractionParticlesPressures.at(i) );
159 velocityDifference += ( currentInteractionParticlesVelocities.at(2 * ( i - 1 ) + 1) - previousInteractionParticlesVelocities.at(2 * ( i - 1 ) + 1) ) * ( currentInteractionParticlesVelocities.at(2 * ( i - 1 ) + 1) - previousInteractionParticlesVelocities.at(2 * ( i - 1 ) + 1) );
160 velocityDifference += ( currentInteractionParticlesVelocities.at(2 * ( i - 1 ) + 2) - previousInteractionParticlesVelocities.at(2 * ( i - 1 ) + 2) ) * ( currentInteractionParticlesVelocities.at(2 * ( i - 1 ) + 2) - previousInteractionParticlesVelocities.at(2 * ( i - 1 ) + 2) );
161 }
162 pressureDifference = sqrt(pressureDifference);
163 velocityDifference = sqrt(velocityDifference);
164
165 if ( iterationNumber > 0 ) {
166 pressureDifference /= previousInteractionParticlesPressures.computeNorm();
167 velocityDifference /= previousInteractionParticlesVelocities.computeNorm();
168 }
169
170 OOFEM_LOG_RELEVANT("%3d %le %le\n", iterationNumber++, pressureDifference, velocityDifference);
171 }
172 if ( iterationNumber > maxiter ) {
173 OOFEM_ERROR("Maximal fluid-structure interface iteration count exceded");
174 }
175 stepN->incrementStateCounter();
176}
177
178void
179FluidStructureProblem :: preInitializeNextStep()
180{
181 for ( auto &emodel: emodelList ) {
182 emodel->preInitializeNextStep();
183 }
184}
185
186#if 0
187void
188FluidStructureProblem :: postInitializeCurrentStep()
189{
190 for ( int i = 1; i <= nModels; i++ ) {
191 DIIDynamic *dynamicProblem = dynamic_cast< DIIDynamic * >( this->giveSlaveProblem(i) );
192 if ( dynamicProblem ) {
193 this->giveCurrentStep()->setTimeDiscretization( dynamicProblem->giveInitialTimeDiscretization() );
194 }
195 }
196}
197#endif
198} // end namespace oofem
#define REGISTER_EngngModel(class)
TimeDiscretizationType giveInitialTimeDiscretization()
Definition diidynamic.h:126
Dof * giveDofWithID(int dofID) const
Definition dofmanager.C:127
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Definition domain.h:461
DofManager * giveDofManager(int n)
Definition domain.C:317
int ndomains
Number of receiver domains.
Definition engngm.h:215
Domain * giveDomain(int n)
Definition engngm.C:1936
double computeNorm() const
Definition floatarray.C:861
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
int maxiter
Max number of iterations.
double rtolv
Convergence tolerance.
int iterationNumber
Iteration counter.
void giveCoupledVelocities(FloatArray &answer, TimeStep *stepN)
TimeDiscretizationType giveInitialTimeDiscretization()
double giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof) override
Definition pfem.C:257
StaggeredProblem(int i, EngngModel *_master=nullptr)
EngngModel * giveSlaveProblem(int i) override
Returns i-th slave problem.
TimeStep * giveCurrentStep(bool force=false) override
void setRenumberFlag() override
Sets the renumber flag to true.
std ::vector< std ::unique_ptr< EngngModel > > emodelList
List of engineering models to solve sequentially.
void incrementStateCounter()
Updates solution state counter.
Definition timestep.h:213
double giveTargetTime()
Returns target time.
Definition timestep.h:164
int giveNumber()
Returns receiver's number.
Definition timestep.h:144
bool isTheFirstStep()
Definition timestep.C:148
#define OOFEM_ERROR(...)
Definition error.h:79
#define _IFT_FluidStructureProblem_maxiter
#define _IFT_FluidStructureProblem_rtolp
#define _IFT_FluidStructureProblem_rtolv
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define OOFEM_LOG_RELEVANT(...)
Definition logger.h:142

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