LifeV
FastAssemblerMixed.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 the LifeV library
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
21  License along with this library; if not, see <http://www.gnu.org/licenses/>
22 
23 
24 *******************************************************************************
25 */
26 //@HEADER
27 
28 /*!
29  * @file
30  @brief This file contains the definition of the FastAssemblerMixed class.
31 
32  This function is used to perform an efficient assembly of FE matrices
33  where the test and trial functions are not the same.
34 
35  @date 06/2016
36  @author Davide Forti <davide.forti@epfl.ch>
37  */
38 
39 #ifndef FASTASSEMBLERMIXED_HPP
40 #define FASTASSEMBLERMIXED_HPP
41 
42 #include <lifev/core/LifeV.hpp>
43 #include <lifev/core/mesh/RegionMesh.hpp>
44 #include <lifev/core/array/VectorEpetra.hpp>
45 #include <lifev/core/array/MatrixEpetra.hpp>
46 #include <lifev/core/fem/QuadratureRule.hpp>
47 #include <lifev/core/fem/ReferenceFE.hpp>
48 #include <lifev/core/fem/CurrentFE.hpp>
49 #include <lifev/core/fem/FESpace.hpp>
50 
51 namespace LifeV
52 {
53 
55 {
56 public:
57 
59  typedef boost::shared_ptr<mesh_Type> meshPtr_Type;
60 
62  typedef boost::shared_ptr<vector_Type> vectorPtr_Type;
63 
65  typedef boost::shared_ptr<matrix_Type> matrixPtr_Type;
66 
68  typedef boost::shared_ptr< comm_Type > commPtr_Type;
69 
71  typedef boost::shared_ptr< qr_Type > qrPtr_Type;
72 
73  typedef FESpace<mesh_Type, MapEpetra> fespace_Type;
74  typedef boost::shared_ptr<fespace_Type> fespacePtr_Type;
75 
76 
77  //! Constructor
78  /*!
79  * @param mesh - input mesh
80  * @param comm - communicator
81  * @param refFE_test - reference FE space of test functions
82  * @param refFE_trial - reference FE space of trial functions
83  * @param qr_integration - quadrature rule to be used for the integration
84  */
85  FastAssemblerMixed( const meshPtr_Type& mesh, const commPtr_Type& comm,
86  const ReferenceFE* refFE_test, const ReferenceFE* refFE_trial,
87  const qr_Type* qr_integration );
88 
89  //! Destructor
91 
92  //! @name Methods
93  //@{
94 
95  //! Allocate space for members before the assembly
96  /*!
97  * @param numElements - data file
98  * @param fe_test - current FE test functions
99  * @param fespace_test - FE space test functions
100  * @param fe_trial - current FE trial functions
101  * @param fespace_trial - FE space trial functions
102  */
103  void allocateSpace( const int& numElements,
104  CurrentFE* fe_test, const fespacePtr_Type& fespace_test,
105  CurrentFE* fe_trial, const fespacePtr_Type& fespace_trial );
106 
107  //! FE Assembly block (0,1) of Navier-Stokes
108  /*!
109  * @param matrix - global matrix, in this case the block (0,1) of Navier-Stokes
110  */
111  void assemble_NS_block01 ( matrixPtr_Type& matrix );
112 
113  //! FE Assembly block (1,0) of Navier-Stokes
114  /*!
115  * @param matrix - global matrix, in this case the block (1,0) of Navier-Stokes
116  */
117  void assemble_NS_block10 ( matrixPtr_Type& matrix );
118 
119  //! FE Assembly block (1,0) of SUPG stabilization for Navier-Stokes
120  /*!
121  * @param matrix - global matrix, in this case the block (1,0) of Navier-Stokes
122  * @param u_h - vector extrapolapolated velocity
123  */
124  void assemble_SUPG_block10 ( matrixPtr_Type& matrix, const vector_Type& u_h );
125 
126  //! FE Assembly block (0,1) of SUPG stabilization for Navier-Stokes
127  /*!
128  * @param matrix - global matrix, in this case the block (0,1) of Navier-Stokes
129  * @param u_h - vector extrapolapolated velocity
130  */
131  void assemble_SUPG_block01 ( matrixPtr_Type& matrix, const vector_Type& u_h );
132 
133  //@}
134 
135 private:
136 
139 
143 
144  double * M_detJacobian;
145  double *** M_invJacobian;
146 
147  double ** M_phi_test;
148  double *** M_dphi_test;
149  double ** M_phi_trial;
150  double *** M_dphi_trial;
151  double ** M_elements_test;
152  double ** M_elements_trial;
153 
157 
158  double**** M_vals;
159  int** M_rows;
160  int** M_cols;
161 
162 };
163 
164 } // Namespace LifeV
165 
166 #endif // FASTASSEMBLERMIXED_HPP
VectorEpetra - The Epetra Vector format Wrapper.
boost::shared_ptr< mesh_Type > meshPtr_Type
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
const ReferenceFE * M_referenceFE_trial
boost::shared_ptr< matrix_Type > matrixPtr_Type
void assemble_SUPG_block01(matrixPtr_Type &matrix, const vector_Type &u_h)
FE Assembly block (0,1) of SUPG stabilization for Navier-Stokes.
void updateInverseJacobian(const UInt &iQuadPt)
RegionMesh< LinearTetra > mesh_Type
void assemble_NS_block10(matrixPtr_Type &matrix)
FE Assembly block (1,0) of Navier-Stokes.
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
void allocateSpace(const int &numElements, CurrentFE *fe_test, const fespacePtr_Type &fespace_test, CurrentFE *fe_trial, const fespacePtr_Type &fespace_trial)
Allocate space for members before the assembly.
void assemble_SUPG_block10(matrixPtr_Type &matrix, const vector_Type &u_h)
FE Assembly block (1,0) of SUPG stabilization for Navier-Stokes.
double Real
Generic real data.
Definition: LifeV.hpp:175
void assemble_NS_block01(matrixPtr_Type &matrix)
FE Assembly block (0,1) of Navier-Stokes.
CurrentFE - A primordial class for the assembly of the local matrices/vectors retaining the values on...
Definition: CurrentFE.hpp:243
boost::shared_ptr< comm_Type > commPtr_Type
The class for a reference Lagrangian finite element.
boost::shared_ptr< qr_Type > qrPtr_Type
const ReferenceFE * M_referenceFE_test
FastAssemblerMixed(const meshPtr_Type &mesh, const commPtr_Type &comm, const ReferenceFE *refFE_test, const ReferenceFE *refFE_trial, const qr_Type *qr_integration)
Constructor.
FESpace< mesh_Type, MapEpetra > fespace_Type
QuadratureRule - The basis class for storing and accessing quadrature rules.
boost::shared_ptr< vector_Type > vectorPtr_Type
boost::shared_ptr< fespace_Type > fespacePtr_Type
MatrixEpetra< Real > matrix_Type