LifeV
FSIcouplingCE.cpp
Go to the documentation of this file.
1 #include <lifev/fsi_blocks/solver/FSIcouplingCE.hpp>
2 
3 namespace LifeV
4 {
5 
6 FSIcouplingCE::FSIcouplingCE ( const commPtr_Type& communicator ) :
7 M_comm ( communicator )
8 {
9 }
10 
12 {
13 }
14 
15 void
16 FSIcouplingCE::setUp ( const Real& timeStep, const Real& interfaceDofs, const Real& beta, const Real& gamma,
17  const mapPtr_Type& interfaceMap, const FESpacePtr_Type& fluidVelocityFESpace,
18  const FESpacePtr_Type& structureDisplacementFESpace, const vectorPtr_Type& numerationInterface)
19 {
20  M_timeStep = timeStep;
21  M_interface = interfaceDofs; // Number of interface dofs per component
22  M_fluidVelocityFESpace = fluidVelocityFESpace;
23  M_structureDisplacementFESpace = structureDisplacementFESpace;
24  M_interfaceMap = interfaceMap;
25  M_numerationInterface = numerationInterface;
26  M_beta = beta;
27  M_gamma = gamma;
28 }
29 
30 void
31 FSIcouplingCE::setUp ( const Real& timeStep, const Real& interfaceDofs, const Real& coefficientBDF,
32  const mapPtr_Type& interfaceMap, const FESpacePtr_Type& fluidVelocityFESpace,
33  const FESpacePtr_Type& structureDisplacementFESpace, const vectorPtr_Type& numerationInterface )
34 {
35  M_timeStep = timeStep;
36  M_interface = interfaceDofs; // Number of interface dofs per component
37  M_fluidVelocityFESpace = fluidVelocityFESpace;
38  M_structureDisplacementFESpace = structureDisplacementFESpace;
39  M_interfaceMap = interfaceMap;
40  M_numerationInterface = numerationInterface;
41  M_coefficientBDF = coefficientBDF;
42 }
43 
44 void
45 FSIcouplingCE::buildBlocks ( std::map<ID, ID> const& localDofMap, const bool& lambda_num_structure, bool useBDF )
46 {
47  // if lambda_num_structure = true, lambda follows the numeration of the structure
48  // if lambda_num_structure = false, lambda follows the numeration of the fluid
49 
50  Real value = 1.0;
51  std::map<ID, ID>::const_iterator ITrow;
52 
53  M_lambdaToFluidMomentum.reset ( new MatrixEpetra<Real> ( M_fluidVelocityFESpace->map() ) );
54 
55  if ( lambda_num_structure )
56  {
57  for (UInt dim = 0; dim < 3; ++dim)
58  {
59  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
60  {
61  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->second) ) >= 0 )
62  {
63  M_lambdaToFluidMomentum->addToCoefficient( ITrow->first + dim * M_fluidVelocityFESpace->dof().numTotalDof(),
64  (int) (*M_numerationInterface)[ITrow->second] + dim * M_interface,
65  value );
66  }
67  }
68  }
69  }
70  else
71  {
72  for (UInt dim = 0; dim < 3; ++dim)
73  {
74  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
75  {
76  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->first) ) >= 0 )
77  {
78  M_lambdaToFluidMomentum->addToCoefficient( ITrow->first + dim * M_fluidVelocityFESpace->dof().numTotalDof(),
79  (int) (*M_numerationInterface)[ITrow->first] + dim * M_interface,
80  value );
81  }
82  }
83  }
84  }
85 
86  M_lambdaToFluidMomentum->globalAssemble(M_interfaceMap, M_fluidVelocityFESpace->mapPtr());
87 
88  // ----------------------------------------------------------------------------------------------------------------------------
89 
90  value = -1.0;
91  M_lambdaToStructureMomentum.reset ( new MatrixEpetra<Real> ( M_structureDisplacementFESpace->map() ) );
92 
93  if ( lambda_num_structure )
94  {
95  for (UInt dim = 0; dim < 3; ++dim)
96  {
97  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
98  {
99  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->second) ) >= 0 )
100  {
101  M_lambdaToStructureMomentum->addToCoefficient( ITrow->second + dim * M_structureDisplacementFESpace->dof().numTotalDof(),
102  (int) (*M_numerationInterface)[ITrow->second] + dim * M_interface,
103  value );
104  }
105  }
106  }
107  }
108  else
109  {
110  for (UInt dim = 0; dim < 3; ++dim)
111  {
112  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
113  {
114  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->first) ) >= 0 )
115  {
116  M_lambdaToStructureMomentum->addToCoefficient( ITrow->second + dim * M_structureDisplacementFESpace->dof().numTotalDof(),
117  (int) (*M_numerationInterface)[ITrow->first] + dim * M_interface,
118  value );
119  }
120  }
121  }
122  }
123 
124  M_lambdaToStructureMomentum->globalAssemble(M_interfaceMap, M_structureDisplacementFESpace->mapPtr());
125 
126  // ----------------------------------------------------------------------------------------------------------------------------
127 
128  value = 1.0;
129  M_fluidVelocityToLambda.reset ( new MatrixEpetra<Real> ( *M_interfaceMap ) );
130 
131  if ( lambda_num_structure )
132  {
133  for (UInt dim = 0; dim < 3; ++dim)
134  {
135  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
136  {
137  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->second) ) >= 0 )
138  {
139  M_fluidVelocityToLambda->addToCoefficient( (int) (*M_numerationInterface)[ITrow->second] + dim * M_interface,
140  ITrow->first + dim * M_fluidVelocityFESpace->dof().numTotalDof(),
141  value );
142  }
143  }
144  }
145  }
146  else
147  {
148  for (UInt dim = 0; dim < 3; ++dim)
149  {
150  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
151  {
152  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->first) ) >= 0 )
153  {
154  M_fluidVelocityToLambda->addToCoefficient( (int) (*M_numerationInterface)[ITrow->first] + dim * M_interface,
155  ITrow->first + dim * M_fluidVelocityFESpace->dof().numTotalDof(),
156  value );
157  }
158  }
159  }
160  }
161 
162  M_fluidVelocityToLambda->globalAssemble(M_fluidVelocityFESpace->mapPtr(), M_interfaceMap);
163 
164  // ----------------------------------------------------------------------------------------------------------------------------
165 
166  if ( useBDF )
167  {
168  value = ( (-1.0) / ( M_timeStep ) * M_coefficientBDF );
169  }
170  else
171  {
172  value = ( (-1.0 * M_gamma) / ( M_timeStep * M_beta ) );
173  }
174 
175  M_structureDisplacementToLambda.reset ( new MatrixEpetra<Real> ( *M_interfaceMap ) );
176 
177  if ( lambda_num_structure )
178  {
179  for (UInt dim = 0; dim < 3; ++dim)
180  {
181  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
182  {
183  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->second) ) >= 0 )
184  {
185  M_structureDisplacementToLambda->addToCoefficient( (int) (*M_numerationInterface)[ITrow->second] + dim * M_interface,
186  ITrow->second + dim * M_structureDisplacementFESpace->dof().numTotalDof(),
187  value );
188  }
189  }
190  }
191  }
192  else
193  {
194  for (UInt dim = 0; dim < 3; ++dim)
195  {
196  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
197  {
198  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->first) ) >= 0 )
199  {
200  M_structureDisplacementToLambda->addToCoefficient( (int) (*M_numerationInterface)[ITrow->first] + dim * M_interface,
201  ITrow->second + dim * M_structureDisplacementFESpace->dof().numTotalDof(),
202  value );
203  }
204  }
205  }
206  }
207 
208  M_structureDisplacementToLambda->globalAssemble(M_structureDisplacementFESpace->mapPtr(), M_interfaceMap);
209 
210  // ----------------------------------------------------------------------------------------------------------------------------
211 
212  value = -1.0;
213  M_structureDisplacementToFluidDisplacement.reset( new MatrixEpetra<Real> ( M_fluidVelocityFESpace->map() ) );
214 
215  if ( lambda_num_structure )
216  {
217  for (UInt dim = 0; dim < 3; ++dim)
218  {
219  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
220  {
221  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->second) ) >= 0 )
222  {
223  M_structureDisplacementToFluidDisplacement->addToCoefficient( ITrow->first + dim * M_fluidVelocityFESpace->dof().numTotalDof(),
224  ITrow->second + dim * M_structureDisplacementFESpace->dof().numTotalDof(),
225  value );
226  }
227  }
228  }
229  }
230  else
231  {
232  for (UInt dim = 0; dim < 3; ++dim)
233  {
234  for (ITrow = localDofMap.begin(); ITrow != localDofMap.end(); ++ITrow)
235  {
236  if ( M_numerationInterface->map().map (Unique)->LID ( static_cast<EpetraInt_Type> (ITrow->first) ) >= 0 )
237  {
238  M_structureDisplacementToFluidDisplacement->addToCoefficient( ITrow->first + dim * M_fluidVelocityFESpace->dof().numTotalDof(),
239  ITrow->second + dim * M_structureDisplacementFESpace->dof().numTotalDof(),
240  value );
241  }
242  }
243  }
244  }
245 
246  M_structureDisplacementToFluidDisplacement->globalAssemble(M_structureDisplacementFESpace->mapPtr(), M_fluidVelocityFESpace->mapPtr());
247 }
248 
249 
250 
251 }
void buildBlocks(std::map< ID, ID > const &locDofMap, const bool &lambda_num_structure, bool useBDF=false)
Builds the coupling blocks.
void updateInverseJacobian(const UInt &iQuadPt)
commPtr_Type M_comm
communicator
FESpacePtr_Type M_structureDisplacementFESpace
std::shared_ptr< VectorEpetra > vectorPtr_Type
~FSIcouplingCE()
Destructor.
void setUp(const Real &timeStep, const Real &interfaceDofs, const Real &coefficientBDF, const mapPtr_Type &interfaceMap, const FESpacePtr_Type &fluidVelocityFESpace, const FESpacePtr_Type &structureDisplacementFESpace, const vectorPtr_Type &numerationInterface)
Set parameters. To be used when Newmark is used on the structure.
FESpacePtr_Type M_fluidVelocityFESpace
std::shared_ptr< map_Type > mapPtr_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
vectorPtr_Type M_numerationInterface
FSIcouplingCE - File handling the coupling blocks when conforming discretizations are used...
void setUp(const Real &timeStep, const Real &interfaceDofs, const Real &beta, const Real &gamma, const mapPtr_Type &interfaceMap, const FESpacePtr_Type &fluidVelocityFESpace, const FESpacePtr_Type &structureDisplacementFESpace, const vectorPtr_Type &numerationInterface)
Set parameters. To be used when Newmark is used on the structure.
std::shared_ptr< comm_Type > commPtr_Type
FSIcouplingCE(const commPtr_Type &communicator)
Constructor.
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
mapPtr_Type M_interfaceMap
std::shared_ptr< FESpace_Type > FESpacePtr_Type