YaDICs  V04.14.a
Yet another Digital Image Correlation software: platform dedicated to 2/3D Fluid and Solid kinematics field measurements.
 All Classes Files Functions Variables Pages
Public Member Functions | Public Attributes | List of all members
Ccorrelation_opticalFlow< T, Timg > Class Template Referenceabstract

This is the mother class of Optical Flow strategies: It implements different steps required for the iterative and parallel resolution of the linear system K.U = F. Notice that the parallelization is done on elements for OFFEM and on pixels for OFI. The function "stiffness" implements the calcul of K, i.e. the global stiffness tensor. It is calculated only one time and since it is based on the gradient of the img[0]. The function "forces" implements the calcul of F, i.e. the global forces vector. It is updated at each loop from img[1]* = img[1](x+U,y+V,z+W). Functions Fassembly and Kassembly are specific to OFFEM or OFI, thus are virtual pure, and implements within deeper class how K anf F assembly are done. The function "gradient" implements the point to point gradient. This function allows specifying how the gradient should be done in different case, i.e. centered, left or right one. The function "actu" implements the field actualization from new U vector as well as the difference between img[0] and img[1]* required for the further calcul of F. The function "regularize" implements a median filter used to enforce the smoothing of the solution vector "U" at each loop. It ensure a better convergency of the iterative calcul and is a first step in the implementation of a penality strategy (to do). The function "solve_cimg" is specific to OFFEM and OFI. Specifically when one use mask. It also allow using different solvers. More...

#include <Ccorrelation_opticalFlow.h>

Inheritance diagram for Ccorrelation_opticalFlow< T, Timg >:
Ccorrelation< T, Timg > Ccorrelation_opticalFlow_fem< T, Timg > Ccorrelation_opticalFlow_fem_gradient< T, Timg > Ccorrelation_opticalFlow_fem_newton< T, Timg > Ccorrelation_opticalFlow_integrated< T, Timg > Ccorrelation_opticalFlow_integrated_block< T, Timg > Ccorrelation_opticalFlow_integrated_block_A< T, Timg > Ccorrelation_opticalFlow_integrated_block_A< T, Timg >

Public Member Functions

virtual void init (CParameterNetCDF &fp)
 
virtual void exec (const Cimage< Timg > &oImage, const Cmesh< T, Timg > &oMesh, Cfield< T, Timg > &oField)
 
virtual void assignG (const Cmesh< T, Timg > &oMesh, CImgList< T > &K, CImgList< T > &Ke)
 
virtual void assignL (const Cmesh< T, Timg > &oMesh, CImgList< T > &Ve)=0
 
virtual void reset (const int &thread, CImgList< T > &Ke)
 
virtual void evalGN (const std::vector< T > &grad, const CImgList< T > &N, CImgList< T > &Ve)=0
 
virtual void e_ref (const Cimage< Timg > &oImage, const Cmesh< T, Timg > &oMesh, const int &elem, std::vector< T > &ref)=0
 
virtual void evalU (const CImgList< T > &solution, const Cmesh< T, Timg > &oMesh, const int &elem, const CImgList< T > &N, std::vector< T > &U)=0
 
virtual void evalShape (const CImg< Timg > &list, const Cimage< Timg > &oImage, const Cfield< T, Timg > &oField, const std::vector< T > &ref, const std::vector< T > &lCoor, const std::vector< T > &gCoor, CImgList< T > &N)=0
 
virtual void gradient (const Cimage< Timg > &oImage, const std::vector< T > &x, std::vector< T > &grad)
 
void evalRes_loop (const Cimage< Timg > &oImage, const std::vector< T > &U, const std::vector< T > &x, T &residue, T &deformed)
 
void evalRes_final (const Cimage< Timg > &oImage, const std::vector< T > &U, const std::vector< T > &x, T &residue, T &deformed)
 
virtual void tensors (const int &thread, const CImgList< T > &Ve, const T &residue, CImgList< T > &Ke)
 
virtual void assembly (const int &thread, const Cmesh< T, Timg > &oMesh, const int &elem, const CImgList< T > &Ke, CImgList< T > &K)=0
 
virtual void metric (const Cimage< Timg > &oImage, const Cmesh< T, Timg > &oMesh, CImgList< T > &solution, Cfield< T, Timg > &oField, CImgList< T > &N, CImgList< T > &Ve, CImgList< T > &Ke, std::vector< T > &ref, std::vector< T > &offset, std::vector< T > &lCoor, std::vector< T > &gCoor, std::vector< T > &U, std::vector< T > &grad, T &residue, T &deformed, CImg< Timg > &list, CImg< int > &p, CImg< T > &current_res, CImgList< T > &K)
 
virtual void solve (CImgList< T > &K, const Cmesh< T, Timg > &oMesh, CImgList< T > &solution)
 
virtual void test (const CImgList< T > &solution, const double &epsilon, Cfield< T, Timg > &oField, T &previous_res, int &iter, int &converged)
 
virtual void export_res (const Cimage< Timg > &oImage, const Cmesh< T, Timg > &oMesh, Cfield< T, Timg > &oField, CImgList< T > &N, std::vector< T > &ref, std::vector< T > &offset, std::vector< T > &lCoor, std::vector< T > &gCoor, std::vector< T > &U, T &residue, T &deformed, CImg< Timg > &list)
 
virtual void store (Cfield< T, Timg > &oField, const Cimage< Timg > &oImage, const std::vector< T > gCoor, const std::vector< T > &U, const T &deformed, const T &residue)
 
- Public Member Functions inherited from Ccorrelation< T, Timg >
virtual void regularize (const Cmesh< T, Timg > &oMesh, const Cimage< Timg > &oImage, CImgList< T > &solution)
 
void count (const int &increment, const int &dim, int &counter)
 

Public Attributes

std::string m_shape_type
 
std::string m_interpLoop
 
std::string m_interpFinal
 
double m_epsilon
 
std::vector< T > m_ref
 
- Public Attributes inherited from Ccorrelation< T, Timg >
std::string class_name
 
std::string m_correl_type
 
std::string m_correl_name
 
m_current_res
 
int m_threadNB
 
int _3D
 
int m_medianFilter
 
bool m_verbose
 

Detailed Description

template<typename T, typename Timg>
class Ccorrelation_opticalFlow< T, Timg >

This is the mother class of Optical Flow strategies: It implements different steps required for the iterative and parallel resolution of the linear system K.U = F. Notice that the parallelization is done on elements for OFFEM and on pixels for OFI. The function "stiffness" implements the calcul of K, i.e. the global stiffness tensor. It is calculated only one time and since it is based on the gradient of the img[0]. The function "forces" implements the calcul of F, i.e. the global forces vector. It is updated at each loop from img[1]* = img[1](x+U,y+V,z+W). Functions Fassembly and Kassembly are specific to OFFEM or OFI, thus are virtual pure, and implements within deeper class how K anf F assembly are done. The function "gradient" implements the point to point gradient. This function allows specifying how the gradient should be done in different case, i.e. centered, left or right one. The function "actu" implements the field actualization from new U vector as well as the difference between img[0] and img[1]* required for the further calcul of F. The function "regularize" implements a median filter used to enforce the smoothing of the solution vector "U" at each loop. It ensure a better convergency of the iterative calcul and is a first step in the implementation of a penality strategy (to do). The function "solve_cimg" is specific to OFFEM and OFI. Specifically when one use mask. It also allow using different solvers.

Parameters
[in]m_shape_type: <int>. This is the shape function type. For OFFEM, Q4 is implemented, and for OFI, rigid_body 2D and 3D, homogeneous 2D and 3D, dilatation 2D and 3D, brasilian 2D, blade 2D are implemented.
[in]m_epsilon: <double>. This is the limit rate of convergence. Under this value calcul is considered as converged within a minimum.

Member Function Documentation

template<typename T , typename Timg >
void Ccorrelation_opticalFlow< T, Timg >::evalRes_final ( const Cimage< Timg > &  oImage,
const std::vector< T > &  U,
const std::vector< T > &  x,
T &  residue,
T &  deformed 
)
inline

the displacement U[dim] in physical unit is divided by the Ratio Pix/millimeter considering the current size of super pixel.//

template<typename T , typename Timg >
void Ccorrelation_opticalFlow< T, Timg >::evalRes_loop ( const Cimage< Timg > &  oImage,
const std::vector< T > &  U,
const std::vector< T > &  x,
T &  residue,
T &  deformed 
)
inline

the displacement U[dim] in physical unit is divided by the Ratio Pix/millimeter considering the current size of super pixel.//

template<typename T , typename Timg >
virtual void Ccorrelation_opticalFlow< T, Timg >::gradient ( const Cimage< Timg > &  oImage,
const std::vector< T > &  x,
std::vector< T > &  grad 
)
inlinevirtual

implements border conditions (global centered gradient, but right and left one for left and right borders respectively). An other check is done for pixels, in an element, but at the border of a mask.

Reimplemented in Ccorrelation_opticalFlow_integrated< T, Timg >.

template<typename T , typename Timg >
virtual void Ccorrelation_opticalFlow< T, Timg >::init ( CParameterNetCDF &  fp)
inlinevirtual

The documentation for this class was generated from the following file: