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
Cfield_opticalFlow_integrated_block.h
1 /**********************************************************************
2  * Copyright (C) 2012, The YaDICs Project Developers.
3  * See the COPYRIGHT file at the top-level directory of this distribution ./COPYRIGHT.
4  * See ./COPYING file for copying and redistribution conditions.
5  *
6  * This file is part of YaDICs.
7  * YaDICs is free software: you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation, either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * YaDICs is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with YaDICs. If not, see <http://www.gnu.org/licenses/>.
19  *
20  * Information about how to use the software are provided at http://yadic.univ-lille1.fr/
21  **********************************************************************/
22 
23 
24 #ifndef CFIELD_OPTICALFLOW_INTEGRATED_BLOCK
25 #define CFIELD_OPTICALFLOW_INTEGRATED_BLOCK
26 
27 
28 //-----------------------------CIMG_NETCDF-------------------------------
29 #define cimg_plugin "plugins/add_fileformat.h"
30 #define cimg_use_netcdf
31 
32 #ifdef cimg_use_netcdf
33  #include "../NetCDF.Tool/struct_parameter_NetCDF.h"
34 // #define cimg_plugin1 "plugins/netcdf_file_format4CImg1.h"
35 // #define cimglist_plugin1 "plugins/netcdf_file_format4CImgList1.h"
36  #define cimg_plugin2 "plugins/netcdf_file_format4CImg2.h"
37  #define cimglist_plugin2 "plugins/netcdf_file_format4CImgList2.h"
38  #include "../CImg.Tool/CImg_NetCDF.h"
39 #endif
40 //-----------------------------CIMG_NETCDF-------------------------------
41 
42 #include "Cfield_opticalFlow_integrated.h"
43 
48 template<typename T, typename Timg>
50 {
51 
52  public:
53 
55  {
57  this->class_name = "Field : Optical Flow Integrated per Block";
58 
59  }
60 
61  virtual void init(CParameterNetCDF &fp)
62  {
65 
66  }
67 
68  virtual void exec(const Cimage<Timg> &oImage, const Cmesh<T,Timg> &oMesh, CImgList<T> &modes_prev, CshapeFunction<T,Timg> &oShape_prev, const std::vector<T> &ref_prev, const CImgList<T> &nodes_prev)
69  {
72  this->m_mode.assign(this->m_pShape->m_dof, oMesh.grid_dims(0),oMesh.grid_dims(1),oMesh.grid_dims(2),1, 0);
73 
74  Cfield<T,Timg>::exec(oImage, oMesh, modes_prev, oShape_prev, ref_prev, nodes_prev);
75 
76  }
77 
78  virtual void project(const Cmesh<T,Timg> &oMesh, const Cimage<Timg> &oImage, CImgList<T> &modes_prev, CshapeFunction<T,Timg> &oShape_prev, const std::vector<T> &ref_prev, const CImgList<T> &nodes_prev, CImgList<T> &modes_cur)
79  {//project modes on new base
80 
81  //Check if rigid body translations are required within the current shape function. If yes, [U,V,W] at nodes are transfered to [Tx,Ty,Tz] at element centroid.
82  std::vector<int> md, pos;
83  md.clear();pos.clear();
84  for (int i = 0; i<this->m_pShape->m_mods.size(); i++)
85  {
86  if (this->m_pShape->m_mods[i].compare(0,2,"Tx") == 0)
87  {pos.push_back(i);md.push_back(0);}
88  if (this->m_pShape->m_mods[i].compare(0,2,"Ty") == 0)
89  {pos.push_back(i);md.push_back(1);}
90  if (this->m_pShape->m_mods[i].compare(0,2,"Tz") == 0)
91  {pos.push_back(i);md.push_back(2);}
92  }
93 
94  if ((oShape_prev.m_name).compare(0,5,"OFFEM") == 0)
95  {//previous step was a FEM type
96 
97  for (int i=0; i<pos.size(); i++)
98  {//linear interpolation of required rigid body motions
99 
100  modes_prev[md[i]].resize(this->grid_dims(0),this->grid_dims(1),this->grid_dims(2),this->grid_dims(3),3);
101  modes_cur[pos[i]] = modes_prev[md[i]];
102 
103  }
104 
105  }
106  else if((oShape_prev.m_name).compare(0,2,"IC") == 0)
107  {
108 
109  //build U map on previous mesh
110  CImgList<T> U(oMesh.num_var(), modes_prev[0].width(), modes_prev[0].height(), modes_prev[0].depth(), 1, 0);
111  cimglist_for(U, mode)
112  {U[mode] = modes_prev[mode+2];}
113 
114  for (int i=0; i<pos.size(); i++)
115  {//linear interpolation of required rigid body motions
116 
117  U[md[i]].resize(this->grid_dims(0),this->grid_dims(1),this->grid_dims(2),this->grid_dims(3),3);
118  modes_cur[pos[i]] = U[md[i]];
119 
120  }
121 
122  }
123  else if((oShape_prev.m_name).compare(0,4,"OFIB") == 0)
124  {
125 
126  //One assumes that modes doesn't change from one OFIB to another
127  if (modes_cur.size() == modes_prev.size())
128  {this->Field2Field(modes_prev, modes_cur);}
129  else
130  {printf("warning : modes have been changes for multi-scale OFIB strategy!!!");}
131 
132  }
133  else if((oShape_prev.m_name).compare(0,3,"OFI") == 0)
134  {
135 
136  //build U map on previous mesh
137  CImgList<T> U(oMesh.num_var(), modes_cur[0].width(), modes_cur[0].height(), modes_cur[0].depth(), 1, 0);
138  this->OFIU(oMesh, oImage, oShape_prev, ref_prev, nodes_prev, modes_prev, U);
139 
140  for (int i=0; i<pos.size(); i++)
141  {//linear interpolation of required rigid body motions
142  modes_cur[pos[i]] = U[md[i]];
143  }
144 
145  }
146 
147  }
148 
149 
150 };
151 
152 #endif
virtual void exec(const Cimage< Timg > &oImage, const Cmesh< T, Timg > &oMesh, CImgList< T > &modes_prev, CshapeFunction< T, Timg > &oShape_prev, const std::vector< T > &ref_prev, const CImgList< T > &nodes_prev)
Definition: Cfield.h:100
This is the mother class of shape function objets (virtual pure).
Definition: CshapeFunction.h:48
Cfield_opticalFlow_integrated_block()
Definition: Cfield_opticalFlow_integrated_block.h:54
In this class, specific projection field method is implemented. In the case of O.F.F.E.M. the node value must be identified through a projection, so "fromField2Mod" implements a parallelized loop over image pixel and identifies modes from a least square method strategy.
Definition: Cfield_opticalFlow_integrated.h:50
virtual void init(CParameterNetCDF &fp)
Definition: Cfield_opticalFlow_integrated_block.h:61
virtual void init(CParameterNetCDF &fp)
Definition: Cfield_opticalFlow_integrated.h:63
Definition: Cmesh.h:61
In this class, specific projection field method is implemented. In the case of O.F.F.E.M. the node value must be identified through a projection, so "fromField2Mod" implements a parallelized loop over image pixel and identifies modes from a least square method strategy.
Definition: Cfield_opticalFlow_integrated_block.h:49
virtual void exec(const Cimage< Timg > &oImage, const Cmesh< T, Timg > &oMesh, CImgList< T > &modes_prev, CshapeFunction< T, Timg > &oShape_prev, const std::vector< T > &ref_prev, const CImgList< T > &nodes_prev)
Definition: Cfield_opticalFlow_integrated_block.h:68
Definition: Cimage.h:57