LifeV
EvaluationReturnAtQuadraturePoints.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 EvaluationInterpolateValue class.
31 
32  @date 06/2011
33  @author Samuel Quinodoz <samuel.quinodoz@epfl.ch>
34  */
35 
36 #ifndef EVALUATION_RETURN_AT_QUAD_PTS_HPP
37 #define EVALUATION_RETURN_AT_QUAD_PTS_HPP
38 
39 #include <lifev/core/LifeV.hpp>
40 
41 #include <lifev/core/array/VectorEpetra.hpp>
42 #include <lifev/core/array/VectorSmall.hpp>
43 
44 #include <lifev/eta/fem/ETCurrentFE.hpp>
45 #include <lifev/eta/fem/ETCurrentFlag.hpp>
46 #include <lifev/eta/fem/ETFESpace.hpp>
47 #include <lifev/core/fem/QuadratureRule.hpp>
48 
49 #include <lifev/eta/expression/ExpressionReturnAtQuadraturePoints.hpp>
50 
51 namespace LifeV
52 {
53 
54 namespace ExpressionAssembly
55 {
56 
57 //! Evaluation for the interpolation of a FE function
58 /*!
59  @author Samuel Quinodoz <samuel.quinodoz@epfl.ch>
60 
61  This class aims at representing an interpolated FE value.
62 
63  This is the generic implementation, so representing a vectorial FE
64 
65  This class is an Evaluation class, and therefore, has all the methods
66  required to work within the Evaluation trees.
67  */
68 template<typename MeshType, typename MapType, UInt SpaceDim, UInt FieldDim>
69 class EvaluationReturnAtQuadraturePoints
70 {
71 public:
72 
73  //! @name Public Types
74  //@{
75 
76  //! Type of the value returned by this class
77  typedef VectorSmall<FieldDim> return_Type;
78 
79  //! Type of the FESpace that has to be used with this class
80  typedef ETFESpace<MeshType, MapType, SpaceDim, FieldDim> fespace_Type;
81 
82  //! Pointer on the FESpace
83  typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
84 
85  //! Vector of the values
86  typedef std::vector<std::vector<return_Type>> vector_Type;
87 
88  //@}
89 
90 
91  //! @name Static constants
92  //@{
93 
94  //! Flag for the global current FE
95  const static flag_Type S_globalUpdateFlag;
96 
97  //! Flag for the test current FE
98  const static flag_Type S_testUpdateFlag;
99 
100  //! Flag for the solution current FE
101  const static flag_Type S_solutionUpdateFlag;
102 
103  //@}
104 
105 
106  //! @name Constructors, destructor
107  //@{
108 
109  //! Copy constructor
110  EvaluationReturnAtQuadraturePoints (const EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, FieldDim>& evaluation)
111  :
112  M_fespace ( evaluation.M_fespace),
113  M_vector ( evaluation.M_vector ),
114  M_quadrature (0),
115  //M_currentFE(M_fespace->refFE(),M_fespace->geoMap()),
116  M_currentFE (evaluation.M_currentFE),
117  M_returnedValues (evaluation.M_returnedValues)
118  {
119  if (evaluation.M_quadrature != 0)
120  {
121  M_quadrature = new QuadratureRule (* (evaluation.M_quadrature) );
122  //M_currentFE.setQuadratureRule(*(evaluation.M_quadrature));
123  }
124  }
125 
126  //! Expression-based constructor
127  explicit EvaluationReturnAtQuadraturePoints (const ExpressionReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, FieldDim>& expression)
128  :
129  M_fespace ( expression.fespace() ),
130  M_vector ( expression.vector() ),
131  M_quadrature (0),
132  M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
133  M_returnedValues (0)
134  {}
135 
136  //! Destructor
137  ~EvaluationReturnAtQuadraturePoints()
138  {
139  if (M_quadrature != 0)
140  {
141  delete M_quadrature;
142  }
143  }
144 
145  //@}
146 
147 
148  //! @name Methods
149  //@{
150 
151  //! Interal update, computes the interpolated values.
152  void update (const UInt& iElement)
153  {
154  zero();
155 
156  for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
157  {
158  M_returnedValues[q] = M_vector[iElement][q];
159  }
160  }
161 
162  //! Erase all the interpolated values stored internally
163  void zero()
164  {
165  for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
166  {
167  for (UInt iDim (0); iDim < FieldDim; ++iDim)
168  {
169  M_returnedValues[q][iDim] = 0.0;
170  }
171  }
172  }
173 
174  //! Show the values interpolated
175  void showValues() const
176  {
177  std::cout << " Values : " << std::endl;
178 
179  for (UInt i (0); i < M_quadrature->nbQuadPt(); ++i)
180  {
181  std::cout << M_returnedValues[i] << std::endl;
182  }
183  }
184 
185  //! Display method
186  static void display (std::ostream& out = std::cout)
187  {
188  out << "returned[" << FieldDim << "]";
189  }
190 
191  //@}
192 
193 
194  //! @name Set Methods
195  //@{
196 
197  //! Do nothing setter for the global current FE
198  template< typename CFEType >
199  void setGlobalCFE (const CFEType* /*globalCFE*/)
200  {}
201 
202  //! Do nothing setter for the test current FE
203  template< typename CFEType >
204  void setTestCFE (const CFEType* /*testCFE*/)
205  {}
206 
207  //! Do nothing setter for the solution current FE
208  template< typename CFEType >
209  void setSolutionCFE (const CFEType* /*solutionCFE*/)
210  {}
211 
212  //! Setter for the quadrature rule (deep copy)
213  void setQuadrature (const QuadratureRule& qr)
214  {
215  if (M_quadrature != 0)
216  {
217  delete M_quadrature;
218  }
219  M_quadrature = new QuadratureRule (qr);
220  M_currentFE.setQuadratureRule (qr);
221  M_returnedValues.resize (qr.nbQuadPt() );
222  }
223 
224  //@}
225 
226 
227  //! @name Get Methods
228  //@{
229 
230  //! Getter for a value
231  return_Type value_q (const UInt& q) const
232  {
233  return M_returnedValues[q];
234  }
235 
236  //! Getter for the value for a vector
237  return_Type value_qi (const UInt& q, const UInt& /*i*/) const
238  {
239  return M_returnedValues[q];
240  }
241 
242  //! Getter for the value for a matrix
243  return_Type value_qij (const UInt& q, const UInt& /*i*/, const UInt& /*j*/) const
244  {
245  return M_returnedValues[q];
246  }
247 
248  //@}
249 
250 private:
251 
252  //! @name Private Methods
253  //@{
254 
255  //! No empty constructor
256  EvaluationReturnAtQuadraturePoints();
257 
258  //@}
259 
260  fespacePtr_Type M_fespace;
261  vector_Type M_vector;
262 
263  QuadratureRule* M_quadrature;
264  ETCurrentFE<SpaceDim, 1> M_currentFE;
265 
266  std::vector<return_Type> M_returnedValues;
267 };
268 
269 
270 //! Evaluation for the interpolation of a FE function
271 /*!
272  @author Samuel Quinodoz <samuel.quinodoz@epfl.ch>
273 
274  This class aims at representing an interpolated FE value.
275 
276  This is the specialized (partially) implementation representing a scalar FE
277 
278  This class is an Evaluation class, and therefore, has all the methods
279  required to work within the Evaluation trees.
280  */
281 template<typename MeshType, typename MapType, UInt SpaceDim>
282 class EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, 1>
283 {
284 
285 public:
286 
287  //! @name Public Types
288  //@{
289 
290  //! Type of the value returned by this class
291  typedef Real return_Type;
292 
293  //! Type of the FESpace to be used in this class
294  typedef ETFESpace<MeshType, MapType, SpaceDim, 1> fespace_Type;
295 
296  //! Type of the pointer on the FESpace
297  typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
298 
299  //! Type of the vector to be used
300  typedef std::vector<std::vector<VectorSmall<1>>> vector_Type;
301 
302  //@}
303 
304 
305  //! @name Static constants
306  //@{
307 
308  //! Flag for the global current FE
309  const static flag_Type S_globalUpdateFlag;
310 
311  //! Flag for the test current FE
312  const static flag_Type S_testUpdateFlag;
313 
314  //! Flag for the solution current FE
315  const static flag_Type S_solutionUpdateFlag;
316 
317  //@}
318 
319 
320  //! @name Constructors, destructor
321  //@{
322 
323  //! Copy constructor
324  EvaluationReturnAtQuadraturePoints (const EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, 1>& evaluation)
325  :
326  M_fespace ( evaluation.M_fespace),
327  M_vector ( evaluation.M_vector ),
328  M_quadrature (0),
329  //M_currentFE(M_fespace->refFE(),M_fespace->geoMap()),
330  M_currentFE (evaluation.M_currentFE),
331  M_returnedValues (evaluation.M_returnedValues)
332  {
333  if (evaluation.M_quadrature != 0)
334  {
335  M_quadrature = new QuadratureRule (* (evaluation.M_quadrature) );
336  //M_currentFE.setQuadratureRule(*(evaluation.M_quadrature));
337  }
338  }
339 
340  //! Expression-based constructor
341  explicit EvaluationReturnAtQuadraturePoints (const ExpressionReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, 1>& expression)
342  :
343  M_fespace ( expression.fespace() ),
344  M_vector ( expression.vector() ),
345  M_quadrature (0),
346  M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
347  M_returnedValues (0)
348  {}
349 
350  //! Destructor
351  ~EvaluationReturnAtQuadraturePoints()
352  {
353  if (M_quadrature != 0)
354  {
355  delete M_quadrature;
356  }
357  }
358 
359  //@}
360 
361 
362  //! @name Methods
363  //@{
364 
365  //! Internal update: computes the interpolated value
366  void update (const UInt& iElement)
367  {
368  zero();
369  for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
370  {
371  M_returnedValues[q] = M_vector[iElement][q](0);
372  }
373  }
374 
375  //! Erase the interpolated values stored internally
376  void zero()
377  {
378  for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
379  {
380  M_returnedValues[q] = 0.0;
381  }
382  }
383 
384  //! Show the values
385  void showValues() const
386  {
387  std::cout << " Values : " << std::endl;
388 
389  for (UInt i (0); i < M_quadrature->nbQuadPt(); ++i)
390  {
391  std::cout << M_returnedValues[i] << std::endl;
392  }
393  }
394 
395  //! Display method
396  static void display (std::ostream& out = std::cout)
397  {
398  out << "returned[1]";
399  }
400 
401  //@}
402 
403 
404  //! @name Set Methods
405  //@{
406 
407  //! Do nothing setter for the global current FE
408  template< typename CFEType >
409  void setGlobalCFE (const CFEType* /*globalCFE*/)
410  {}
411 
412  //! Do nothing setter for the test current FE
413  template< typename CFEType >
414  void setTestCFE (const CFEType* /*testCFE*/)
415  {}
416 
417  //! Do nothing setter for the solution current FE
418  template< typename CFEType >
419  void setSolutionCFE (const CFEType* /*solutionCFE*/)
420  {}
421 
422  //! Setter for the quadrature rule
423  void setQuadrature (const QuadratureRule& qr)
424  {
425  if (M_quadrature != 0)
426  {
427  delete M_quadrature;
428  }
429  M_quadrature = new QuadratureRule (qr);
430  M_currentFE.setQuadratureRule (qr);
431  M_returnedValues.resize (qr.nbQuadPt() );
432  }
433 
434  //@}
435 
436 
437  //! @name Get Methods
438  //@{
439 
440  //! Getter for a value
441  return_Type value_q (const UInt& q) const
442  {
443  return M_returnedValues[q];
444  }
445 
446  //! Getter for the value for a vector
447  return_Type value_qi (const UInt& q, const UInt& /*i*/) const
448  {
449  return M_returnedValues[q];
450  }
451 
452  //! Getter for the value for a matrix
453  return_Type value_qij (const UInt& q, const UInt& /*i*/, const UInt& /*j*/) const
454  {
455  return M_returnedValues[q];
456  }
457 
458  //@}
459 
460 private:
461 
462  //! @name Private Methods
463  //@{
464 
465  //! No empty constructor
466  EvaluationReturnAtQuadraturePoints();
467 
468  //@}
469 
470  //! Data storage
471  fespacePtr_Type M_fespace;
472  vector_Type M_vector;
473  QuadratureRule* M_quadrature;
474 
475  //! Structure for the computations
476  ETCurrentFE<SpaceDim, 1> M_currentFE;
477 
478  //! Storage for the temporary values
479  std::vector<return_Type> M_returnedValues;
480 };
481 
482 
483 
484 template<typename MeshType, typename MapType, UInt SpaceDim, UInt FieldDim>
485 const flag_Type
486 EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, FieldDim>::
487 S_globalUpdateFlag = ET_UPDATE_NONE;
488 
489 template<typename MeshType, typename MapType, UInt SpaceDim, UInt FieldDim>
490 const flag_Type
491 EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, FieldDim>::
492 S_testUpdateFlag = ET_UPDATE_NONE;
493 
494 template<typename MeshType, typename MapType, UInt SpaceDim, UInt FieldDim>
495 const flag_Type
496 EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, FieldDim>::
497 S_solutionUpdateFlag = ET_UPDATE_NONE;
498 
499 
500 
501 template<typename MeshType, typename MapType, UInt SpaceDim>
502 const flag_Type
503 EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, 1>::
504 S_globalUpdateFlag = ET_UPDATE_NONE;
505 
506 template<typename MeshType, typename MapType, UInt SpaceDim>
507 const flag_Type
508 EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, 1>::
509 S_testUpdateFlag = ET_UPDATE_NONE;
510 
511 template<typename MeshType, typename MapType, UInt SpaceDim>
512 const flag_Type
513 EvaluationReturnAtQuadraturePoints<MeshType, MapType, SpaceDim, 1>::
514 S_solutionUpdateFlag = ET_UPDATE_NONE;
515 
516 } // Namespace ExpressionAssembly
517 
518 } // Namespace LifeV
519 #endif
uint32_type flag_Type
bit-flag with up to 32 different flags
Definition: LifeV.hpp:197
QuadratureRule(const QuadratureRule &qr)
Copy constructor.
void updateInverseJacobian(const UInt &iQuadPt)
double Real
Generic real data.
Definition: LifeV.hpp:175
QuadratureRule - The basis class for storing and accessing quadrature rules.
const UInt & nbQuadPt() const
Getter for the number of quadrature points.
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191