LifeV
MultiscaleModelMultiscale.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 Model Multiscale
30  *
31  * @date 12-03-2009
32  * @author Cristiano Malossi <cristiano.malossi@epfl.ch>
33  *
34  * @maintainer Cristiano Malossi <cristiano.malossi@epfl.ch>
35  */
36 
37 #ifndef MultiscaleModelMultiscale_H
38 #define MultiscaleModelMultiscale_H 1
39 
40 #include <lifev/multiscale/framework/MultiscaleCommunicatorsManager.hpp>
41 
42 #include <lifev/multiscale/algorithms/MultiscaleAlgorithm.hpp>
43 
44 #include <lifev/multiscale/couplings/MultiscaleCoupling.hpp>
45 
46 #include <lifev/multiscale/models/MultiscaleModel.hpp>
47 
48 #if defined(LIFEV_HAS_ZERODIMENSIONAL)
49 #include <lifev/multiscale/models/MultiscaleModelWindkessel0D.hpp>
50 #include <lifev/multiscale/models/MultiscaleModel0D.hpp>
51 #endif
52 
53 #if defined(LIFEV_HAS_ONEDFSI)
54 #include <lifev/multiscale/models/MultiscaleModelFSI1D.hpp>
55 #endif
56 
57 #if defined(LIFEV_HAS_NAVIERSTOKES)
58 #include <lifev/multiscale/models/MultiscaleModelFluid3D.hpp>
59 #endif
60 
61 #if defined(LIFEV_HAS_FSI)
62 #include <lifev/multiscale/models/MultiscaleModelFSI3D.hpp>
63 #endif
64 
65 namespace LifeV
66 {
67 namespace Multiscale
68 {
69 
70 //! MultiscaleModelMultiscale - Multiscale model
71 /*!
72  * @author Cristiano Malossi
73  *
74  * @see Full description of the Geometrical Multiscale Framework: \cite Malossi-Thesis
75  * @see Methodology: \cite Malossi2011Algorithms \cite Malossi2011Algorithms1D \cite Malossi2011Algorithms3D1DFSI \cite BlancoMalossi2012
76  * @see Applications: \cite Malossi2011Algorithms3D1DFSIAortaIliac \cite LassilaMalossi2012IdealLeftVentricle \cite BonnemainMalossi2012LVAD
77  *
78  * The MultiscaleModelMultiscale class is an implementation of the MultiscaleModel
79  * for a general multiscale problem.
80  */
82 {
83 public:
84 
85  //! @name Constructors & Destructor
86  //@{
87 
88  //! Constructor
89  explicit MultiscaleModelMultiscale();
90 
91  //! Destructor
92  virtual ~MultiscaleModelMultiscale();
93 
94  //@}
95 
96 
97  //! @name MultiscaleModel Methods
98  //@{
99 
100  //! Setup the data of the model.
101  /*!
102  * @param fileName Name of data file.
103  */
104  void setupData ( const std::string& fileName );
105 
106  //! Setup the model.
107  void setupModel();
108 
109  //! Build the initial model.
110  void buildModel();
111 
112  //! Update the model.
113  void updateModel();
114 
115  //! Solve the model.
116  void solveModel();
117 
118  //! Update the solution.
119  void updateSolution();
120 
121  //! Save the solution
122  void saveSolution();
123 
124  //! Display some information about the multiscale problem (call after SetupProblem)
125  void showMe();
126 
127  //! Return a specific scalar quantity to be used for a comparison with a reference value.
128  /*!
129  * This method is meant to be used for night checks.
130  * @return reference quantity.
131  */
132  Real checkSolution() const;
133 
134  //@}
135 
136 
137  //! @name Methods
138  //@{
139 
140  //! Build the global map for the coupling vectors
141  /*!
142  * @param couplingMap Global coupling map
143  */
144  void createCouplingMap ( MapEpetra& couplingMap );
145 
146  //! Import the values of the coupling variables
147  void importCouplingVariables ( const multiscaleVector_Type& couplingVariables );
148 
149  //! Export the values of the coupling variables
150  void exportCouplingVariables ( multiscaleVector_Type& couplingVariables );
151 
152  //! Compute the values of the interface residuals
154 
155  //! Export the values of the interface residuals
156  void exportCouplingResiduals ( multiscaleVector_Type& couplingResiduals );
157 
158  //! Export the Jacobian matrix
159  /*!
160  * @param jacobian Matrix
161  */
162  void exportJacobian ( multiscaleMatrix_Type& jacobian );
163 
164  //! Check if the topology is changed
165  /*!
166  * A topology change can be caused by a change in the coupling equations by,
167  * for example, the opening/closure of a valve (see MultiscaleCouplingFlowRateValve).
168  * @return true if the topology is changed, false otherwise
169  */
170  bool topologyChange();
171 
172  //@}
173 
174 
175  //! @name Get Methods
176  //@{
177 
178  //! Get the number of the coupling variables
179  /*!
180  * @return number of the coupling variables
181  */
183 
184  //@}
185 
186 private:
187 
188  //! @name Unimplemented Methods
189  //@{
190 
192 
194 
195  //@}
196 
197  // Models & Couplings
199 
200  // Models & Couplings
203 
204  // Algorithm for subiterations
206 };
207 
208 //! Factory create function
210 {
211  return new MultiscaleModelMultiscale();
212 }
213 
214 } // Namespace multiscale
215 } // Namespace LifeV
216 
217 #endif /* MultiscaleModelMultiscale_H */
void exportCouplingVariables(multiscaleVector_Type &couplingVariables)
Export the values of the coupling variables.
std::shared_ptr< multiscaleAlgorithm_Type > multiscaleAlgorithmPtr_Type
void showMe()
Display some information about the multiscale problem (call after SetupProblem)
bool topologyChange()
Check if the topology is changed.
Real checkSolution() const
Return a specific scalar quantity to be used for a comparison with a reference value.
MultiscaleModelMultiscale(const MultiscaleModelMultiscale &model)
void updateInverseJacobian(const UInt &iQuadPt)
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
MultiscaleModel multiscaleModel_Type
void createCouplingMap(MapEpetra &couplingMap)
Build the global map for the coupling vectors.
MatrixEpetra< Real > multiscaleMatrix_Type
std::vector< multiscaleModelPtr_Type > multiscaleModelsContainer_Type
std::vector< multiscaleCouplingPtr_Type > multiscaleCouplingsContainer_Type
void computeCouplingResiduals()
Compute the values of the interface residuals.
double Real
Generic real data.
Definition: LifeV.hpp:175
void exportJacobian(multiscaleMatrix_Type &jacobian)
Export the Jacobian matrix.
void exportCouplingResiduals(multiscaleVector_Type &couplingResiduals)
Export the values of the interface residuals.
MultiscaleModelMultiscale - Multiscale model.
UInt couplingVariablesNumber()
Get the number of the coupling variables.
void importCouplingVariables(const multiscaleVector_Type &couplingVariables)
Import the values of the coupling variables.
multiscaleModel_Type * createMultiscaleModelMultiscale()
Factory create function.
MultiscaleCommunicatorsManager - The Multiscale Communicators Manager.
MultiscaleModelMultiscale & operator=(const MultiscaleModelMultiscale &model)
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
VectorEpetra multiscaleVector_Type