36 #ifndef EVALUATION_INTERPOLATE_VALUE_HPP 37 #define EVALUATION_INTERPOLATE_VALUE_HPP 39 #include <lifev/core/LifeV.hpp> 41 #include <lifev/core/array/VectorEpetra.hpp> 42 #include <lifev/core/array/VectorSmall.hpp> 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> 49 #include <lifev/eta/expression/ExpressionInterpolateValue.hpp> 51 #include <boost/shared_ptr.hpp> 57 namespace ExpressionAssembly
71 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
72 class EvaluationInterpolateValue
83 typedef ETFESpace<MeshType, MapType, SpaceDim, FieldDim> fespace_Type;
86 typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
113 EvaluationInterpolateValue (
const EvaluationInterpolateValue<MeshType, MapType, SpaceDim, FieldDim>& evaluation)
115 M_fespace ( evaluation.M_fespace),
116 M_vector ( evaluation.M_vector, Repeated),
119 M_currentFE (evaluation.M_currentFE),
120 M_interpolatedValues (evaluation.M_interpolatedValues)
122 if (evaluation.M_quadrature != 0)
130 explicit EvaluationInterpolateValue (
const ExpressionInterpolateValue<MeshType, MapType, SpaceDim, FieldDim>& expression)
132 M_fespace ( expression.fespace() ),
133 M_vector ( expression.vector(), Repeated ),
135 M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
136 M_interpolatedValues (0)
140 ~EvaluationInterpolateValue()
142 if (M_quadrature != 0)
155 void update (
const UInt& iElement)
159 for (UInt i (0); i < M_fespace->refFE().nbDof(); ++i)
161 for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
163 for (UInt iDim (0); iDim < FieldDim; ++iDim)
165 UInt globalID (M_fespace->dof().localToGlobalMap (iElement, i)
166 + iDim * M_fespace->dof().numTotalDof() );
168 M_interpolatedValues[q][iDim] +=
169 M_currentFE.phi (i, q)
170 * M_vector[globalID];
181 for (
UInt iDim (0); iDim < FieldDim; ++iDim)
183 M_interpolatedValues[q][iDim] = 0.0;
189 void showValues()
const 191 std::cout <<
" Values : " << std::endl;
195 std::cout << M_interpolatedValues[i] << std::endl;
200 static void display (std::ostream& out = std::cout)
202 out <<
"interpolated[" << FieldDim <<
"]";
212 template<
typename CFEType >
213 void setGlobalCFE (
const CFEType* )
217 template<
typename CFEType >
218 void setTestCFE (
const CFEType* )
222 template<
typename CFEType >
223 void setSolutionCFE (
const CFEType* )
229 if (M_quadrature != 0)
234 M_currentFE.setQuadratureRule (qr);
235 M_interpolatedValues.resize (qr.nbQuadPt() );
245 return_Type value_q (
const UInt& q)
const 247 return M_interpolatedValues[q];
251 return_Type value_qi (
const UInt& q,
const UInt& )
const 253 return M_interpolatedValues[q];
257 return_Type value_qij (
const UInt& q,
const UInt& ,
const UInt& )
const 259 return M_interpolatedValues[q];
270 EvaluationInterpolateValue();
274 fespacePtr_Type M_fespace;
275 vector_Type M_vector;
278 ETCurrentFE<SpaceDim, 1> M_currentFE;
280 std::vector<return_Type> M_interpolatedValues;
295 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
296 class EvaluationInterpolateValue<MeshType, MapType, SpaceDim, 1>
308 typedef ETFESpace<MeshType, MapType, SpaceDim, 1> fespace_Type;
311 typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
338 EvaluationInterpolateValue (
const EvaluationInterpolateValue<MeshType, MapType, SpaceDim, 1>& evaluation)
340 M_fespace ( evaluation.M_fespace),
341 M_vector ( evaluation.M_vector, Repeated),
344 M_currentFE (evaluation.M_currentFE),
345 M_interpolatedValues (evaluation.M_interpolatedValues)
347 if (evaluation.M_quadrature != 0)
355 explicit EvaluationInterpolateValue (
const ExpressionInterpolateValue<MeshType, MapType, SpaceDim, 1>& expression)
357 M_fespace ( expression.fespace() ),
358 M_vector ( expression.vector(), Repeated ),
360 M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
361 M_interpolatedValues (0)
365 ~EvaluationInterpolateValue()
367 if (M_quadrature != 0)
380 void update (
const UInt& iElement)
384 for (UInt i (0); i < M_fespace->refFE().nbDof(); ++i)
386 for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
388 UInt globalID (M_fespace->dof().localToGlobalMap (iElement, i) );
390 M_interpolatedValues[q] +=
391 M_currentFE.phi (i, q)
392 * M_vector[globalID];
402 M_interpolatedValues[q] = 0.0;
407 void showValues()
const 409 std::cout <<
" Values : " << std::endl;
413 std::cout << M_interpolatedValues[i] << std::endl;
418 static void display (std::ostream& out = std::cout)
420 out <<
"interpolated[1]";
430 template<
typename CFEType >
431 void setGlobalCFE (
const CFEType* )
435 template<
typename CFEType >
436 void setTestCFE (
const CFEType* )
440 template<
typename CFEType >
441 void setSolutionCFE (
const CFEType* )
447 if (M_quadrature != 0)
452 M_currentFE.setQuadratureRule (qr);
453 M_interpolatedValues.resize (qr.nbQuadPt() );
463 return_Type value_q (
const UInt& q)
const 465 return M_interpolatedValues[q];
469 return_Type value_qi (
const UInt& q,
const UInt& )
const 471 return M_interpolatedValues[q];
475 return_Type value_qij (
const UInt& q,
const UInt& ,
const UInt& )
const 477 return M_interpolatedValues[q];
488 EvaluationInterpolateValue();
493 fespacePtr_Type M_fespace;
494 vector_Type M_vector;
498 ETCurrentFE<SpaceDim, 1> M_currentFE;
501 std::vector<return_Type> M_interpolatedValues;
506 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
508 EvaluationInterpolateValue<MeshType, MapType, SpaceDim, FieldDim>::
509 S_globalUpdateFlag = ET_UPDATE_NONE;
511 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
513 EvaluationInterpolateValue<MeshType, MapType, SpaceDim, FieldDim>::
514 S_testUpdateFlag = ET_UPDATE_NONE;
516 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
518 EvaluationInterpolateValue<MeshType, MapType, SpaceDim, FieldDim>::
519 S_solutionUpdateFlag = ET_UPDATE_NONE;
523 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
525 EvaluationInterpolateValue<MeshType, MapType, SpaceDim, 1>::
526 S_globalUpdateFlag = ET_UPDATE_NONE;
528 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
530 EvaluationInterpolateValue<MeshType, MapType, SpaceDim, 1>::
531 S_testUpdateFlag = ET_UPDATE_NONE;
533 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
535 EvaluationInterpolateValue<MeshType, MapType, SpaceDim, 1>::
536 S_solutionUpdateFlag = ET_UPDATE_NONE;
VectorEpetra - The Epetra Vector format Wrapper.
uint32_type flag_Type
bit-flag with up to 32 different flags
QuadratureRule(const QuadratureRule &qr)
Copy constructor.
void updateInverseJacobian(const UInt &iQuadPt)
double Real
Generic real data.
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)