LifeV
MultiscaleCouplingMeanTotalNormalStressArea.cpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  * @file
29  * @brief File containing the multiscale mean normal stress coupling class
30  *
31  * @date 11-10-2012
32  * @author Cristiano Malossi <cristiano.malossi@epfl.ch>
33  *
34  * @maintainer Cristiano Malossi <cristiano.malossi@epfl.ch>
35  */
36 
37 #include <lifev/multiscale/couplings/MultiscaleCouplingMeanTotalNormalStressArea.hpp>
38 
39 namespace LifeV
40 {
41 namespace Multiscale
42 {
43 
44 // ===================================================
45 // Constructors & Destructor
46 // ===================================================
49 {
50 
51 #ifdef HAVE_LIFEV_DEBUG
52  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::MultiscaleCouplingMeanTotalNormalStressArea() \n";
53 #endif
54 
56 }
57 
58 // ===================================================
59 // Multiscale PhysicalCoupling Implementation
60 // ===================================================
61 void
63 {
64 
65 #ifdef HAVE_LIFEV_DEBUG
66  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::setupCouplingVariablesNumber() \n";
67 #endif
68 
70 }
71 
72 void
74 {
75 
76 #ifdef HAVE_LIFEV_DEBUG
77  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::setupCoupling() \n";
78 #endif
79 
80  super_Type::setupCoupling();
81 
82  // Preliminary checks
83  if ( myModelsNumber() > 0 )
84  {
85  if ( modelsNumber() > 2 )
86  {
87  std::cerr << "!!! ERROR: MultiscaleCouplingMeanTotalNormalStressArea does not work with more than two models !!!" << std::endl;
88  }
89 
90  //TODO: add a check on the type of the two models: they must be a FSI3D and a 1D model.
91  }
92 
93  if ( myModelsNumber() > 0 )
94  {
95  // Impose area boundary condition on the FSI3D model
96  for ( UInt i ( 0 ); i < 2; ++i )
97  if ( myModel ( i ) )
98  if ( M_models[i]->type() == FSI3D )
99  {
100  M_localCouplingFunctions.push_back ( MultiscaleCouplingFunction ( this, 3 ) );
101  multiscaleDynamicCast< MultiscaleInterface > ( M_models[i] )->imposeBoundaryArea ( M_boundaryIDs[i], std::bind ( &MultiscaleCouplingFunction::function, M_localCouplingFunctions.back(), std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, std::placeholders::_4, std::placeholders::_5 ) );
102 
103  break;
104  }
105  }
106 }
107 
108 void
110 {
111 
112 #ifdef HAVE_LIFEV_DEBUG
113  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::initializeCouplingVariables() \n";
114 #endif
115 
116  super_Type::initializeCouplingVariables();
117 
118  // Local variables initialization
119  Real localSum ( 0 );
120  Real globalSum ( 0 );
121 
122  // Compute the area coupling variable as an average of the two areas
123  for ( UInt i ( 0 ); i < 2; ++i )
124  if ( myModel ( i ) )
125  {
126  Real myValue = multiscaleDynamicCast< MultiscaleInterface > ( M_models[i] )->boundaryArea ( M_boundaryIDs[i] );
127  if ( isModelLeaderProcess ( i ) )
128  {
129  localSum += myValue;
130  }
131  }
132 
133  M_comm->SumAll ( &localSum, &globalSum, 1 );
134  if ( myModelsNumber() > 0 )
135  {
136  localCouplingVariables ( 0 ) [3] = globalSum / 2;
137  }
138 }
139 
140 void
142 {
143 
144 #ifdef HAVE_LIFEV_DEBUG
145  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::computeCouplingResiduals() \n";
146 #endif
147 
148  super_Type::computeCouplingResiduals();
149 
150  if ( myModelsNumber() > 0 )
151  {
152  // Impose area boundary condition on the FSI3D model
153  for ( UInt i ( 0 ); i < 2; ++i )
154  if ( myModel ( i ) )
155  if ( M_models[i]->type() == FSI1D )
156  {
157  Real myValueArea = multiscaleDynamicCast< MultiscaleInterface > ( M_models[i] )->boundaryArea ( M_boundaryIDs[i] );
158  if ( isModelLeaderProcess ( i ) )
159  {
160  ( *M_localCouplingResiduals ) [3] = myValueArea - localCouplingVariables ( 0 ) [3];
161  }
162  }
163  }
164 }
165 
166 // ===================================================
167 // Private MultiscaleCoupling Implementation
168 // ===================================================
169 void
171 {
172 
173 #ifdef HAVE_LIFEV_DEBUG
174  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::exportListOfPerturbedModels( localCouplingVariableID ) \n";
175 #endif
176 
177  if ( localCouplingVariableID == 3 )
178  {
179  for ( UInt i ( 0 ); i < 2; ++i )
180  if ( myModel ( i ) )
181  if ( M_models[i]->type() == FSI3D )
182  {
183  perturbedModelsList.push_back ( M_models[i] );
184  }
185  }
186  else
187  {
188  super_Type::exportListOfPerturbedModels ( localCouplingVariableID, perturbedModelsList );
189  }
190 }
191 
192 void
194 {
195 
196 #ifdef HAVE_LIFEV_DEBUG
197  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::insertJacobianConstantCoefficients( jacobian ) \n";
198 #endif
199 
200  super_Type::insertJacobianConstantCoefficients ( jacobian );
201 
202  // The constant coefficients are added by the leader process of model 0.
203  if ( myModel ( 0 ) )
204  if ( isModelLeaderProcess ( 0 ) )
205  {
206  jacobian.addToCoefficient ( M_couplingVariablesOffset + 3, M_couplingVariablesOffset + 3, -1 );
207  }
208 }
209 
210 void
211 MultiscaleCouplingMeanTotalNormalStressArea::insertJacobianDeltaCoefficients ( multiscaleMatrix_Type& jacobian, const UInt& column, const UInt& ID, bool& solveLinearSystem )
212 {
213 
214 #ifdef HAVE_LIFEV_DEBUG
215  debugStream ( 8260 ) << "MultiscaleCouplingMeanTotalNormalStressArea::insertJacobianDeltaCoefficients( jacobian, column, ID, solveLinearSystem ) \n";
216 #endif
217 
218  super_Type::insertJacobianDeltaCoefficients ( jacobian, column, ID, solveLinearSystem );
219 
220  // Model global to local conversion
221  UInt modelLocalID = modelGlobalToLocalID ( ID );
222  if ( myModel ( modelLocalID ) )
223  if ( M_models[modelLocalID]->type() == FSI1D )
224  {
225  Real row ( 0 );
226  Real coefficient ( 0 );
227 
228  row = M_couplingVariablesOffset + 3;
229  coefficient = multiscaleDynamicCast< MultiscaleInterface > ( M_models[modelLocalID] )->boundaryDeltaArea ( M_boundaryIDs[modelLocalID], solveLinearSystem );
230 
231  // Add the coefficient to the matrix
232  if ( isModelLeaderProcess ( modelLocalID ) )
233  {
234  jacobian.addToCoefficient ( row, column, coefficient );
235 
236 #ifdef HAVE_LIFEV_DEBUG
237  debugStream ( 8260 ) << "J(" << row << "," << column << ") = " << coefficient << "\n";
238 #endif
239  }
240  }
241 
242 }
243 
244 } // Namespace Multiscale
245 } // Namespace LifeV
void insertJacobianDeltaCoefficients(multiscaleMatrix_Type &jacobian, const UInt &column, const UInt &ID, bool &solveLinearSystem)
Insert the Jacobian coefficient(s) depending on a perturbation of the model, due to a specific variab...
void insertJacobianConstantCoefficients(multiscaleMatrix_Type &jacobian)
Insert constant coefficients into the Jacobian matrix.
MultiscaleCoupling multiscaleCoupling_Type
void updateInverseJacobian(const UInt &iQuadPt)
MultiscaleCouplingMeanTotalNormalStressArea - Mean normal stress with area coupling condition...
void initializeCouplingVariables()
Initialize the values of the coupling variables.
void exportListOfPerturbedModels(const UInt &localCouplingVariableID, multiscaleModelsContainer_Type &perturbedModelsList)
Build the list of models affected by the perturbation of a local coupling variable.
MatrixEpetra< Real > multiscaleMatrix_Type
std::vector< multiscaleModelPtr_Type > multiscaleModelsContainer_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191