36 #ifndef EVALUATION_INTERPOLATE_GRADIENT_HPP 37 #define EVALUATION_INTERPOLATE_GRADIENT_HPP 39 #include <lifev/core/LifeV.hpp> 41 #include <lifev/core/array/VectorEpetra.hpp> 42 #include <lifev/core/array/VectorSmall.hpp> 43 #include <lifev/core/array/MatrixSmall.hpp> 45 #include <lifev/eta/fem/ETCurrentFE.hpp> 46 #include <lifev/eta/fem/ETCurrentFlag.hpp> 47 #include <lifev/eta/fem/ETFESpace.hpp> 48 #include <lifev/core/fem/QuadratureRule.hpp> 50 #include <lifev/eta/expression/ExpressionInterpolateGradient.hpp> 52 #include <boost/shared_ptr.hpp> 58 namespace ExpressionAssembly
74 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
75 class EvaluationInterpolateGradient
83 typedef MatrixSmall<FieldDim, SpaceDim> return_Type;
86 typedef ETFESpace<MeshType, MapType, SpaceDim, FieldDim> fespace_Type;
89 typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
116 EvaluationInterpolateGradient (
const EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, FieldDim>& evaluation)
118 M_fespace ( evaluation.M_fespace),
119 M_vector ( evaluation.M_vector, Repeated),
121 M_currentFE (evaluation.M_currentFE),
122 M_interpolatedGradients (evaluation.M_interpolatedGradients)
124 if (evaluation.M_quadrature != 0)
131 explicit EvaluationInterpolateGradient (
const ExpressionInterpolateGradient<MeshType, MapType, SpaceDim, FieldDim>& expression)
133 M_fespace ( expression.fespace() ),
134 M_vector ( expression.vector(), Repeated ),
136 M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
137 M_interpolatedGradients (0)
141 ~EvaluationInterpolateGradient()
143 if (M_quadrature != 0)
155 void update (
const UInt& iElement)
159 M_currentFE.update (M_fespace->mesh()->element (iElement), ET_UPDATE_DPHI);
161 for (UInt i (0); i < M_fespace->refFE().nbDof(); ++i)
163 for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
165 for (UInt iDim (0); iDim < SpaceDim; ++iDim)
167 for (UInt jDim (0); jDim < FieldDim; ++jDim)
169 UInt globalID (M_fespace->dof().localToGlobalMap (iElement, i)
170 + jDim * M_fespace->dof().numTotalDof() );
172 M_interpolatedGradients[q][jDim][iDim] +=
173 M_currentFE.dphi (jDim * M_currentFE.nbFEDof() + i, jDim, iDim, q)
174 * M_vector[globalID];
186 for (
UInt i (0); i < SpaceDim; ++i)
188 for (
UInt j (0); j < FieldDim; ++j)
190 M_interpolatedGradients[q][j][i] = 0.0;
197 void showValues()
const 199 std::cout <<
" Gradients : " << std::endl;
203 std::cout << M_interpolatedGradients[i] << std::endl;
208 static void display (std::ostream& out = std::cout)
210 out <<
"interpolated[ " << FieldDim <<
" ][ " << SpaceDim <<
" ]";
220 template<
typename CFEType >
221 void setGlobalCFE (
const CFEType* )
225 template<
typename CFEType >
226 void setTestCFE (
const CFEType* )
230 template<
typename CFEType >
231 void setSolutionCFE (
const CFEType* )
237 if (M_quadrature != 0)
242 M_currentFE.setQuadratureRule (qr);
243 M_interpolatedGradients.resize (qr.nbQuadPt() );
253 return_Type value_q (
const UInt& q)
const 255 return M_interpolatedGradients[q];
259 return_Type value_qi (
const UInt& q,
const UInt& )
const 261 return M_interpolatedGradients[q];
265 return_Type value_qij (
const UInt& q,
const UInt& ,
const UInt& )
const 267 return M_interpolatedGradients[q];
279 EvaluationInterpolateGradient();
284 fespacePtr_Type M_fespace;
285 vector_Type M_vector;
289 ETCurrentFE<SpaceDim, FieldDim> M_currentFE;
292 std::vector<return_Type> M_interpolatedGradients;
295 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
297 EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, FieldDim>::
298 S_globalUpdateFlag = ET_UPDATE_NONE;
300 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
302 EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, FieldDim>::
303 S_testUpdateFlag = ET_UPDATE_NONE;
305 template<
typename MeshType,
typename MapType,
UInt SpaceDim,
UInt FieldDim>
307 EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, FieldDim>::
308 S_solutionUpdateFlag = ET_UPDATE_NONE;
322 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
323 class EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, 1>
335 typedef ETFESpace<MeshType, MapType, SpaceDim, 1> fespace_Type;
338 typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
365 EvaluationInterpolateGradient (
const EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, 1>& evaluation)
367 M_fespace ( evaluation.M_fespace),
368 M_vector ( evaluation.M_vector, Repeated),
369 M_offset ( evaluation.M_offset ),
371 M_currentFE (evaluation.M_currentFE),
372 M_interpolatedGradients (evaluation.M_interpolatedGradients)
374 if (evaluation.M_quadrature != 0)
381 explicit EvaluationInterpolateGradient (
const ExpressionInterpolateGradient<MeshType, MapType, SpaceDim, 1>& expression)
383 M_fespace ( expression.fespace() ),
384 M_vector ( expression.vector(), Repeated ),
385 M_offset ( expression.offset() ),
387 M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
388 M_interpolatedGradients (0)
392 ~EvaluationInterpolateGradient()
394 if (M_quadrature != 0)
407 void update (
const UInt& iElement)
411 M_currentFE.update (M_fespace->mesh()->element (iElement), ET_UPDATE_DPHI);
413 for (UInt i (0); i < M_fespace->refFE().nbDof(); ++i)
415 for (UInt q (0); q < M_quadrature->nbQuadPt(); ++q)
417 UInt globalID (M_fespace->dof().localToGlobalMap (iElement, i) + M_offset);
419 for (UInt iDim (0); iDim < SpaceDim; ++iDim)
421 M_interpolatedGradients[q][iDim] +=
422 M_currentFE.dphi (i, iDim, q)
423 * M_vector[globalID];
434 for (
UInt i (0); i < SpaceDim; ++i)
436 M_interpolatedGradients[q][i] = 0.0;
442 void showValues()
const 444 std::cout <<
" Gradients : " << std::endl;
448 std::cout << M_interpolatedGradients[i] << std::endl;
453 static void display (std::ostream& out = std::cout)
455 out <<
"interpolated[ " << SpaceDim <<
" ]";
465 template<
typename CFEType >
466 void setGlobalCFE (
const CFEType* )
470 template<
typename CFEType >
471 void setTestCFE (
const CFEType* )
475 template<
typename CFEType >
476 void setSolutionCFE (
const CFEType* )
482 if (M_quadrature != 0)
487 M_currentFE.setQuadratureRule (qr);
488 M_interpolatedGradients.resize (qr.nbQuadPt() );
498 return_Type value_q (
const UInt& q)
const 500 return M_interpolatedGradients[q];
504 return_Type value_qi (
const UInt& q,
const UInt& )
const 506 return M_interpolatedGradients[q];
510 return_Type value_qij (
const UInt& q,
const UInt& ,
const UInt& )
const 512 return M_interpolatedGradients[q];
523 EvaluationInterpolateGradient();
528 fespacePtr_Type M_fespace;
529 vector_Type M_vector;
534 ETCurrentFE<SpaceDim, 1> M_currentFE;
537 std::vector<return_Type> M_interpolatedGradients;
541 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
543 EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, 1>::
544 S_globalUpdateFlag = ET_UPDATE_NONE;
546 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
548 EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, 1>::
549 S_testUpdateFlag = ET_UPDATE_NONE;
551 template<
typename MeshType,
typename MapType,
UInt SpaceDim>
553 EvaluationInterpolateGradient<MeshType, MapType, SpaceDim, 1>::
554 S_solutionUpdateFlag = ET_UPDATE_NONE;
570 template<
typename MeshType,
typename MapType>
571 class EvaluationInterpolateGradient<MeshType, MapType, 3, 3>
580 typedef MatrixSmall<3, 3> return_Type;
583 typedef ETFESpace<MeshType, MapType, 3, 3> fespace_Type;
586 typedef std::shared_ptr<fespace_Type> fespacePtr_Type;
613 EvaluationInterpolateGradient (
const EvaluationInterpolateGradient<MeshType, MapType, 3, 3>& evaluation)
615 M_fespace ( evaluation.M_fespace),
616 M_vector ( evaluation.M_vector, Repeated),
617 M_offset ( evaluation.M_offset ),
619 M_currentFE (evaluation.M_currentFE),
620 M_interpolatedGradients (evaluation.M_interpolatedGradients)
622 if (evaluation.M_quadrature != 0)
629 explicit EvaluationInterpolateGradient (
const ExpressionInterpolateGradient<MeshType, MapType, 3, 3>& expression)
631 M_fespace ( expression.fespace() ),
632 M_vector ( expression.vector(), Repeated ),
633 M_offset ( expression.offset() ),
635 M_currentFE (M_fespace->refFE(), M_fespace->geoMap() ),
636 M_interpolatedGradients (0)
640 ~EvaluationInterpolateGradient()
642 if (M_quadrature != 0)
655 void update (
const UInt& iElement)
659 M_currentFE.update (M_fespace->mesh()->element (iElement), ET_UPDATE_DPHI);
660 Real nbFEDof (M_fespace->refFE().nbDof() );
663 MatrixSmall<3, 3> nodalGradMatrix;
665 for (
UInt i (0); i < nbFEDof; ++i)
667 for (
UInt iField (0); iField < 3; ++iField)
669 UInt globalID (M_fespace->dof().localToGlobalMap (iElement, i) + iField * M_fespace->dof().numTotalDof() + M_offset);
670 nodalValues[iField] = M_vector[globalID];
677 nodalGradMatrix = M_currentFE.dphi (i, q)
678 + M_currentFE.dphi (i + nbFEDof, q)
679 + M_currentFE.dphi (i + 2 * nbFEDof, q);
681 M_interpolatedGradients[q] += nodalGradMatrix.emult (nodalValues);
693 for (
UInt j (0); j < 3; ++j)
695 for (
UInt i (0); i < 3; ++i)
697 M_interpolatedGradients[q][j][i] = 0.0;
704 void showValues()
const 706 std::cout <<
" Gradients : " << std::endl;
710 std::cout << M_interpolatedGradients[i] << std::endl;
715 static void display (std::ostream& out = std::cout)
717 out <<
"interpolated[ " << 3 <<
" ]";
727 template<
typename CFEType >
728 void setGlobalCFE (
const CFEType* )
732 template<
typename CFEType >
733 void setTestCFE (
const CFEType* )
737 template<
typename CFEType >
738 void setSolutionCFE (
const CFEType* )
744 if (M_quadrature != 0)
749 M_currentFE.setQuadratureRule (qr);
750 M_interpolatedGradients.resize (qr.nbQuadPt() );
760 return_Type value_q (
const UInt& q)
const 762 return M_interpolatedGradients[q];
766 return_Type value_qi (
const UInt& q,
const UInt& )
const 768 return M_interpolatedGradients[q];
772 return_Type value_qij (
const UInt& q,
const UInt& ,
const UInt& )
const 774 return M_interpolatedGradients[q];
785 EvaluationInterpolateGradient();
790 fespacePtr_Type M_fespace;
791 vector_Type M_vector;
796 ETCurrentFE<3, 3> M_currentFE;
799 std::vector<return_Type> M_interpolatedGradients;
803 template<
typename MeshType,
typename MapType >
805 EvaluationInterpolateGradient<MeshType, MapType, 3, 3>::
806 S_globalUpdateFlag = ET_UPDATE_NONE;
808 template<
typename MeshType,
typename MapType>
810 EvaluationInterpolateGradient<MeshType, MapType, 3, 3>::
811 S_testUpdateFlag = ET_UPDATE_NONE;
813 template<
typename MeshType,
typename MapType>
815 EvaluationInterpolateGradient<MeshType, MapType, 3, 3>::
816 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)