LifeV
MultiscaleCoupling.hpp
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 Physical Coupling
30  *
31  * @date 02-09-2009
32  * @author Cristiano Malossi <cristiano.malossi@epfl.ch>
33  *
34  * @maintainer Cristiano Malossi <cristiano.malossi@epfl.ch>
35  */
36 
37 #ifndef MultiscaleCoupling_H
38 #define MultiscaleCoupling_H 1
39 
40 #include <lifev/multiscale/framework/MultiscaleDefinitions.hpp>
41 #include <lifev/multiscale/framework/MultiscaleGlobalData.hpp>
42 #include <lifev/multiscale/models/MultiscaleModel.hpp>
43 
44 #include <lifev/multiscale/framework/MultiscaleInterface.hpp> // This should not be here
45 
46 namespace LifeV
47 {
48 namespace Multiscale
49 {
50 
51 // Forward declaration
53 
54 //! MultiscaleCoupling - The Multiscale Physical Coupling
55 /*!
56  * @author Cristiano Malossi
57  *
58  * @see Full description of the Geometrical Multiscale Framework: \cite Malossi-Thesis
59  * @see Methodology: \cite Malossi2011Algorithms \cite Malossi2011Algorithms1D \cite Malossi2011Algorithms3D1DFSI \cite BlancoMalossi2012
60  * @see Applications: \cite Malossi2011Algorithms3D1DFSIAortaIliac \cite LassilaMalossi2012IdealLeftVentricle \cite BonnemainMalossi2012LVAD
61  *
62  * The MultiscaleCoupling class provides a general interface between the
63  * MS_Algorithm and all the coupling conditions.
64  */
66 {
67 public:
68 
69  //! @name Type definitions
70  //@{
71 
75 
77 
79 
80  typedef multiscaleVector_Type::combineMode_Type combineMode_Type;
81 
82  //@}
83 
84 
85  //! @name Constructors & Destructor
86  //@{
87 
88  //! Constructor
89  explicit MultiscaleCoupling();
90 
91  //! Destructor
92  virtual ~MultiscaleCoupling() {}
93 
94  //@}
95 
96 
97  //! @name Multiscale PhysicalCoupling Virtual Methods
98  //@{
99 
100  //! Setup the data of the coupling.
101  /*!
102  * @param FileName Name of data file
103  */
104  virtual void setupData ( const std::string& fileName );
105 
106  //! Setup the coupling variables number.
107  virtual void setupCouplingVariablesNumber() = 0;
108 
109  //! Setup the coupling
110  virtual void setupCoupling() = 0;
111 
112  //! Initialize the values of the coupling variables
113  virtual void initializeCouplingVariables() = 0;
114 
115  //! Update the coupling
116  /*!
117  * This method is the analogous of the "updateModel" for the models.
118  * It is alternative to initializeCouplingVariables and is called from the second timestep.
119  * It is reserved for the update of:
120  * <ol>
121  * <li> objects that are not constant with respect to the time but should not be updated during subiterations.
122  * </ol>
123  */
124  virtual void updateCoupling() = 0;
125 
126  //! Compute the values of the local coupling residuals
127  virtual void computeCouplingResiduals() = 0;
128 
129  //! Check if the topology is changed
130  /*!
131  * A topology change can be caused by a change in the coupling equations by,
132  * for example, the opening/closure of a valve (see MultiscaleCouplingFlowRateValve).
133  *
134  * @return true if the topology is changed, false otherwise
135  */
136  virtual bool topologyChange()
137  {
138  return false;
139  }
140 
141  //@}
142 
143 
144  //! @name Methods
145  //@{
146 
147  //! Determine the number of models owned by this coupling
148  /*!
149  * @return number of models owned by the coupling
150  */
151  UInt myModelsNumber() const;
152 
153  //! Determine if the model is owned by this coupling
154  /*!
155  * Note: this method does not check if M_models is empty or not!
156  *
157  * @param localModelID local ID of the model.
158  * @return true if the model is owned by the coupling, false otherwise
159  */
160  bool myModel ( const UInt& localModelID ) const
161  {
162  return M_models[localModelID].get() ? true : false;
163  }
164 
165  //! Determine if this is the model leader process
166  /*!
167  * Note: this method does not check if the model is owned by the process!
168  * Use myModelsNumber() for that!
169  *
170  * @param localModelID local ID of the model.
171  * @return true if this is the model leader process, false otherwise
172  */
173  bool isModelLeaderProcess ( const UInt& localModelID ) const;
174 
175  //! Build the global map for the coupling vectors
176  /*!
177  * @param couplingMap Global coupling map
178  */
179  void createCouplingMap ( MapEpetra& couplingMap );
180 
181  //! Import the values of the coupling variables
182  /*!
183  * @param couplingVariables Global vector of coupling variables
184  */
185  void importCouplingVariables ( const multiscaleVector_Type& couplingVariables )
186  {
187  importCouplingVector ( localCouplingVariables ( 0 ), couplingVariables, Add );
188  }
189 
190  //! Export the values of the coupling variables
191  /*!
192  * @param couplingVariables Global vector of coupling variables
193  */
194  void exportCouplingVariables ( multiscaleVector_Type& couplingVariables )
195  {
196  exportCouplingVector ( couplingVariables, localCouplingVariables ( 0 ), Zero );
197  }
198 
199  //! Export the values of the coupling variables
200  /*!
201  * @param couplingVariables Global vector of coupling variables
202  */
203  void exportCouplingResiduals ( multiscaleVector_Type& couplingResiduals )
204  {
205  exportCouplingVector ( couplingResiduals, *M_localCouplingResiduals, Add );
206  }
207 
208  //! Extrapolate the values of the coupling variables for the next time step
210 
211  //! Lagrange interpolation/extrapolation of the coupling variables at selected time.
212  /*!
213  * @param t interpolation time
214  * @param interpolatedCouplingVariables variables interpolated/extrapolated at time t
215  */
216  void interpolateCouplingVariables ( const Real& t, multiscaleVector_Type& interpolatedCouplingVariables ) const;
217 
218  //! Find if a perturbation is imposed on the coupling.
219  /*!
220  * @return true if a perturbation is imposed
221  */
222  bool isPerturbed() const
223  {
224  return M_perturbedCoupling == -1 ? false : true;
225  }
226 
227  //! Export the Jacobian matrix
228  /*!
229  * @param Jacobian Jacobian Matrix
230  */
231  void exportJacobian ( multiscaleMatrix_Type& jacobian );
232 
233  //! save the coupling variables information on a file
234  void saveSolution();
235 
236  //! Display some information about the coupling
237  void showMe();
238 
239  //! Display the local residuals vector
240  void showMeResiduals() const;
241 
242  //! Display the local coupling variables
243  void showMeCouplingVariables() const;
244 
245  //! Clear the list of pointers to the models.
246  /*!
247  * This method has to be called before the automatic destructor, in order
248  * to disconnect the coupling classes from the model classes.
249  */
251  {
252  M_models.clear();
253  }
254 
255  //@}
256 
257 
258  //! @name Set Methods
259  //@{
260 
261  //! Set the global ID of the coupling condition
262  /*!
263  * @param ID Coupling global ID
264  */
265  void setID ( const UInt& ID )
266  {
267  M_ID = ID;
268  }
269 
270  //! Set the number of models coupled by this coupling condition
271  /*!
272  * @param modelsNumber number of models coupled by this coupling
273  */
274  void setModelsNumber ( const UInt& modelsNumber )
275  {
276  M_models.resize ( modelsNumber );
277  M_boundaryIDs.resize ( modelsNumber );
278  }
279 
280  //! Add a pointer to one of the models to be coupled
281  /*!
282  * @param localModelID local model ID
283  * @param model shared_ptr of the model
284  */
285  void setModel ( const UInt& localModelID, const multiscaleModelPtr_Type& model )
286  {
287  M_models[localModelID] = model;
288  }
289 
290  //! Set the boundary ID of one of the coupled models
291  /*!
292  * @param modelLocalID model local ID
293  * @param boundaryID boundary ID of the model
294  */
295  void setBoundaryID ( const UInt& modelLocalID, const multiscaleID_Type& boundaryLocalID )
296  {
297  M_boundaryIDs[modelLocalID] = boundaryLocalID;
298  }
299 
300  //! Setup the global data of the coupling.
301  /*!
302  * In particular, it can be used to replace the local values specified in
303  * the model data file, with the ones in the global container.
304  *
305  * @param globalData Global data container.
306  */
307  void setGlobalData ( const multiscaleDataPtr_Type& globalData )
308  {
309  M_globalData = globalData;
310  }
311 
312  //! Set the epetra communicator for the coupling
313  /*!
314  * @param comm Epetra communicator
315  */
317  {
318  M_comm = comm;
319  }
320 
321  //@}
322 
323 
324  //! @name Get Methods
325  //@{
326 
327  //! Get the global ID of the coupling
328  /*!
329  * @return global ID of the coupling
330  */
331  const UInt& ID() const
332  {
333  return M_ID;
334  }
335 
336  //! Get the type of the coupling
337  /*!
338  * @return type of the coupling
339  */
340  const couplings_Type& type() const
341  {
342  return M_type;
343  }
344 
345  //! Get the name of the coupling
346  /*!
347  * @return name of the coupling
348  */
349  const std::string& couplingName() const
350  {
351  return M_couplingName;
352  }
353 
354  //! Get the number of models connected by the coupling
355  /*!
356  * @return number of models connected by the coupling
357  */
359  {
360  return M_models.size();
361  }
362 
363  //! Get the model local ID through global ID
364  /*!
365  * @param ID global ID of the model
366  * @return local ID of the model
367  */
368  UInt modelGlobalToLocalID ( const UInt& ID ) const;
369 
370  //! Get the model connected by the coupling through local ID
371  /*!
372  * @param localModelID local ID of the model
373  * @return Pointer to the model
374  */
375  multiscaleModelPtr_Type model ( const UInt& localModelID ) const
376  {
377  return M_models[localModelID];
378  }
379 
380  //! Get the model connected by the coupling through local ID
381  /*!
382  * @param localModelID local ID of the model
383  * @return boundary ID of the model
384  */
385  const multiscaleID_Type& boundaryID ( const UInt& localModelID ) const
386  {
387  return M_boundaryIDs[localModelID];
388  }
389 
390  //! Get the number of the coupling variables
391  /*!
392  * @return number of the coupling variables
393  */
395  {
397  }
398 
399  //! Get the container of the local coupling variables
400  /*!
401  * @return container of the local coupling variables
402  */
404  {
406  }
407 
408  //! Get the perturbed coupling.
409  /*!
410  * If it is unperturbed it returns -1.
411  * @return the localID of the perturbed coupling
412  */
413  const Int& perturbedCoupling() const
414  {
415  return M_perturbedCoupling;
416  }
417 
418  //! Get the local residual.
419  /*!
420  * @return the local residual of the coupling
421  */
423  {
424  return *M_localCouplingResiduals;
425  }
426 
427  //! Get the time interpolation order.
428  /*!
429  * @return the value of the time interpolation order.
430  */
432  {
434  }
435 
436 
437  //@}
438 
439 protected:
440 
441  //! @name Protected MultiscaleCoupling Virtual Methods
442  //@{
443 
444  //! Build the list of models affected by the perturbation of the local coupling variable
445  /*!
446  * @param localCouplingVariableID local coupling variable (perturbed)
447  * @param perturbedModelsList list of models affected by the perturbation of the given coupling variable
448  */
449  virtual void exportListOfPerturbedModels ( const UInt& localCouplingVariableID, multiscaleModelsContainer_Type& perturbedModelsList ) = 0;
450 
451  //! Insert constant coefficients into the Jacobian matrix
452  /*!
453  * @param Jacobian the Jacobian matrix
454  */
455  virtual void insertJacobianConstantCoefficients ( multiscaleMatrix_Type& jacobian ) = 0;
456 
457  //! Insert the Jacobian coefficient(s) depending on a perturbation of the model, due to a specific variable (the column)
458  /*!
459  * @param Jacobian the Jacobian matrix
460  * @param Column the column related to the perturbed variable
461  * @param ID the global ID of the model which is perturbed by the variable
462  * @param SolveLinearSystem a flag to which determine if the linear system has to be solved
463  */
464  virtual void insertJacobianDeltaCoefficients ( multiscaleMatrix_Type& jacobian, const UInt& column, const UInt& ID, bool& linearSystemSolved ) = 0;
465 
466  //@}
467 
468 
469  //! @name Protected Methods
470  //@{
471 
472  //! Access by reference to a specific local coupling variable
473  /*!
474  * This method is used to simplify the access to a specific local coupling variables vector.
475  * Note that the returned value is not const!
476  * @param id id of the local coupling variables vector
477  * @return reference to the local coupling variables vector
478  */
480  {
481  return *M_localCouplingVariables[id];
482  }
484  {
485  return *M_localCouplingVariables[id];
486  }
487 
488  //! Create the local vectors of the coupling
489  void createLocalVectors();
490 
491  //! Reset the history of the couplings
492  /*!
493  * This method is used when the topology change.
494  */
495  void resetCouplingHistory();
496 
497  //! Import the content of the unique global vector into the repeated local vector
498  /*!
499  * @param repeatedLocalVector the repeated local vector
500  * @param uniqueglobalVector the unique global vector
501  */
502  void importCouplingVector ( multiscaleVector_Type& repeatedLocalVector, const multiscaleVector_Type& uniqueGlobalVector, const combineMode_Type& combineMode = Add );
503 
504  //! Export the content of the repeated local vector into the unique global vector
505  /*!
506  * @param uniqueGlobalVector the unique global vector
507  * @param repeatedLocalVector the repeated local vector
508  */
509  void exportCouplingVector ( multiscaleVector_Type& uniqueGlobalVector, const multiscaleVector_Type& repeatedLocalVector, const combineMode_Type& combineMode = Add );
510 
511  //! Display and error message for the specific model
512  /*!
513  * @param model shared_ptr to the specific model
514  */
515  void switchErrorMessage ( const multiscaleModelPtr_Type& model );
516 
517  //@}
518 
521 
523  std::string M_couplingName;
525 
527 
530 
534 
537 
539 
541 
542 private:
543 
544  //! @name Unimplemented Methods
545  //@{
546 
547  MultiscaleCoupling ( const MultiscaleCoupling& coupling );
548 
549  MultiscaleCoupling& operator= ( const MultiscaleCoupling& coupling );
550 
551  //@}
552 };
553 
554 // ===================================================
555 // Inline Methods
556 // ===================================================
557 inline UInt
559 {
560  UInt myModelsNumber ( 0 );
561 
562  for ( UInt i ( 0 ); i < modelsNumber(); ++i )
563  if ( myModel ( i ) )
564  {
565  ++myModelsNumber;
566  }
567 
568  return myModelsNumber;
569 }
570 
571 //! MultiscaleCouplingFunction - The multiscale function for the couplings
572 /*!
573  * @author Cristiano Malossi
574  *
575  * This simple class provides the implementation for the BC function used by the couplings.
576  */
578 {
579 public:
580 
581  //! @name Constructors & Destructor
582  //@{
583 
584  //! Constructor
586 
587  //! Constructor
588  /*!
589  * @param coupling pointer to the coupling
590  * @param id id of the coupling variable
591  */
592  explicit MultiscaleCouplingFunction ( const multiscaleCoupling_Type* coupling, const UInt& id ) : M_coupling ( coupling ), M_couplingID ( id ) {}
593 
594  //! Destructor
596  {
597  /* M_coupling is deleted outside */
598  }
599 
600  //@}
601 
602 
603  //! @name Methods
604  //@{
605 
606  //! Evaluate the coupling quantity
607  /*!
608  * @return evaluation of the function
609  */
610  Real function ( const Real& t, const Real&, const Real&, const Real&, const UInt& )
611  {
612  multiscaleVector_Type interpolatedCouplingVariables ( *M_coupling->couplingVariables() [0] );
613  M_coupling->interpolateCouplingVariables ( t, interpolatedCouplingVariables );
614 
615  return interpolatedCouplingVariables[M_couplingID];
616  }
617 
618  //@}
619 
620 private:
621 
624 };
625 
626 } // Namespace Multiscale
627 } // Namespace LifeV
628 
629 #endif /* MultiscaleCoupling_H */
void setModelsNumber(const UInt &modelsNumber)
Set the number of models coupled by this coupling condition.
void importCouplingVariables(const multiscaleVector_Type &couplingVariables)
Import the values of the coupling variables.
Epetra_CombineMode combineMode_Type
UInt modelGlobalToLocalID(const UInt &ID) const
Get the model local ID through global ID.
virtual void setupCoupling()=0
Setup the coupling.
bool isPerturbed() const
Find if a perturbation is imposed on the coupling.
void switchErrorMessage(const multiscaleModelPtr_Type &model)
Display and error message for the specific model.
UInt myModelsNumber() const
Determine the number of models owned by this coupling.
void setID(const UInt &ID)
Set the global ID of the coupling condition.
MultiscaleCouplingFunction - The multiscale function for the couplings.
void setGlobalData(const multiscaleDataPtr_Type &globalData)
Setup the global data of the coupling.
void showMeCouplingVariables() const
Display the local coupling variables.
Displayer::commPtr_Type multiscaleCommPtr_Type
MultiscaleCouplingFunction couplingFunction_Type
void interpolateCouplingVariables(const Real &t, multiscaleVector_Type &interpolatedCouplingVariables) const
Lagrange interpolation/extrapolation of the coupling variables at selected time.
void createLocalVectors()
Create the local vectors of the coupling.
virtual void updateCoupling()=0
Update the coupling.
std::shared_ptr< multiscaleData_Type > multiscaleDataPtr_Type
void extrapolateCouplingVariables()
Extrapolate the values of the coupling variables for the next time step.
void importCouplingVector(multiscaleVector_Type &repeatedLocalVector, const multiscaleVector_Type &uniqueGlobalVector, const combineMode_Type &combineMode=Add)
Import the content of the unique global vector into the repeated local vector.
void exportCouplingResiduals(multiscaleVector_Type &couplingResiduals)
Export the values of the coupling variables.
MultiscaleCoupling multiscaleCoupling_Type
void setBoundaryID(const UInt &modelLocalID, const multiscaleID_Type &boundaryLocalID)
Set the boundary ID of one of the coupled models.
void exportCouplingVariables(multiscaleVector_Type &couplingVariables)
Export the values of the coupling variables.
const multiscaleID_Type & boundaryID(const UInt &localModelID) const
Get the model connected by the coupling through local ID.
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
void exportCouplingVector(multiscaleVector_Type &uniqueGlobalVector, const multiscaleVector_Type &repeatedLocalVector, const combineMode_Type &combineMode=Add)
Export the content of the repeated local vector into the unique global vector.
const couplingVariablesContainer_Type & couplingVariables() const
Get the container of the local coupling variables.
void updateInverseJacobian(const UInt &iQuadPt)
UInt modelsNumber() const
Get the number of models connected by the coupling.
std::shared_ptr< couplingFunction_Type > couplingFunctionPtr_Type
std::shared_ptr< multiscaleModel_Type > multiscaleModelPtr_Type
virtual bool topologyChange()
Check if the topology is changed.
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
void setCommunicator(const multiscaleCommPtr_Type &comm)
Set the epetra communicator for the coupling.
const UInt & ID() const
Get the global ID of the coupling.
virtual void insertJacobianConstantCoefficients(multiscaleMatrix_Type &jacobian)=0
Insert constant coefficients into the Jacobian matrix.
couplingVariablesContainer_Type M_localCouplingVariables
MultiscaleCouplingFunction(const multiscaleCoupling_Type *coupling, const UInt &id)
Constructor.
void saveSolution()
save the coupling variables information on a file
multiscaleModelsContainer_Type M_models
multiscaleVector_Type & localCouplingVariables(const UInt &id)
std::vector< couplingFunction_Type > couplingFunctionsContainer_Type
const multiscaleVector_Type & residual() const
Get the local residual.
MatrixEpetra< Real > multiscaleMatrix_Type
virtual void setupCouplingVariablesNumber()=0
Setup the coupling variables number.
std::vector< multiscaleModelPtr_Type > multiscaleModelsContainer_Type
bool isModelLeaderProcess(const UInt &localModelID) const
Determine if this is the model leader process.
MultiscaleCoupling(const MultiscaleCoupling &coupling)
void resetCouplingHistory()
Reset the history of the couplings.
const UInt & timeInterpolationOrder() const
Get the time interpolation order.
void exportJacobian(multiscaleMatrix_Type &jacobian)
Export the Jacobian matrix.
couplingFunctionsContainer_Type M_localCouplingFunctions
multiscaleVectorPtr_Type M_localCouplingResiduals
const UInt & couplingVariablesNumber() const
Get the number of the coupling variables.
double Real
Generic real data.
Definition: LifeV.hpp:175
multiscaleIDContainer_Type M_boundaryIDs
void showMeResiduals() const
Display the local residuals vector.
std::vector< multiscaleID_Type > multiscaleIDContainer_Type
virtual void exportListOfPerturbedModels(const UInt &localCouplingVariableID, multiscaleModelsContainer_Type &perturbedModelsList)=0
Build the list of models affected by the perturbation of the local coupling variable.
Real function(const Real &t, const Real &, const Real &, const Real &, const UInt &)
Evaluate the coupling quantity.
const Int & perturbedCoupling() const
Get the perturbed coupling.
multiscaleModelPtr_Type model(const UInt &localModelID) const
Get the model connected by the coupling through local ID.
bool myModel(const UInt &localModelID) const
Determine if the model is owned by this coupling.
void clearModelsList()
Clear the list of pointers to the models.
virtual void insertJacobianDeltaCoefficients(multiscaleMatrix_Type &jacobian, const UInt &column, const UInt &ID, bool &linearSystemSolved)=0
Insert the Jacobian coefficient(s) depending on a perturbation of the model, due to a specific variab...
const std::string & couplingName() const
Get the name of the coupling.
virtual void computeCouplingResiduals()=0
Compute the values of the local coupling residuals.
const couplings_Type & type() const
Get the type of the coupling.
const multiscaleVector_Type & localCouplingVariables(const UInt &id) const
Access by reference to a specific local coupling variable.
void setModel(const UInt &localModelID, const multiscaleModelPtr_Type &model)
Add a pointer to one of the models to be coupled.
void createCouplingMap(MapEpetra &couplingMap)
Build the global map for the coupling vectors.
MultiscaleCoupling - The Multiscale Physical Coupling.
std::vector< multiscaleVectorPtr_Type > couplingVariablesContainer_Type
std::shared_ptr< multiscaleVector_Type > multiscaleVectorPtr_Type
MultiscaleCoupling & operator=(const MultiscaleCoupling &coupling)
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
VectorEpetra multiscaleVector_Type
data_type & operator[](const UInt row)
Access operators.
void showMe()
Display some information about the coupling.
virtual void initializeCouplingVariables()=0
Initialize the values of the coupling variables.