Orfeo Toolbox  3.16
itkFEMSolver.cxx
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: $RCSfile: itkFEMSolver.cxx,v $
5  Language: C++
6  Date: $Date: 2010-03-30 15:20:02 $
7  Version: $Revision: 1.59 $
8 
9  Copyright (c) Insight Software Consortium. All rights reserved.
10  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
11 
12  This software is distributed WITHOUT ANY WARRANTY; without even
13  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
14  PURPOSE. See the above copyright notices for more information.
15 
16 =========================================================================*/
17 
18 // disable debug warnings in MS compiler
19 #ifdef _MSC_VER
20 #pragma warning(disable: 4786)
21 #endif
22 
23 #include "itkFEMSolver.h"
24 
25 #include "itkFEMLoadNode.h"
26 #include "itkFEMLoadElementBase.h"
27 #include "itkFEMElementBase.h"
28 #include "itkFEMLoadBC.h"
29 #include "itkFEMLoadBCMFC.h"
30 #include "itkFEMLoadLandmark.h"
31 
32 #include "itkImageRegionIterator.h"
33 
34 #include <algorithm>
35 
36 namespace itk {
37 namespace fem {
38 
39 /*
40  * Default constructor for Solver class
41  */
42 Solver::Solver() : NGFN(0), NMFC(0)
43 {
45 }
46 
47 void Solver::Clear( void )
48 {
49  this->el.clear();
50  this->node.clear();
51  this->mat.clear();
52  this->load.clear();
53 
54  this->NGFN=0;
55  this->NMFC=0;
57 }
58 
59 
65 {
66  m_ls=ls; // update the pointer to LinearSystemWrapper object
67 
69 }
70 
71 
73 {
74  // set the maximum number of matrices and vectors that
75  // we will need to store inside.
79 }
80 
84 void Solver::Read(std::istream& f) {
85 
86  // clear all arrays - uncomment these 4 lines if you want to use the
87  // Windows-based visualization
88  // el.clear();
89  // node.clear();
90  // mat.clear();
91  // load.clear();
92 
93  // Initialize the pointers to arrays in ReadInfoType object to the
94  // arrays in solver object.
95  ReadInfoType info(&this->node,&this->el,&this->mat);
96 
98  /* then we start reading objects from stream */
99  do
100  {
107  if (!o) { continue; }
108 
113  if ( Node::Pointer o1=dynamic_cast<Node*>(&*o) )
114  {
115  node.push_back(FEMP<Node>(o1));
116  continue;
117  }
118  if ( Material::Pointer o1=dynamic_cast<Material*>(&*o) )
119  {
120  mat.push_back(FEMP<Material>(o1));
121  continue;
122  }
123  if ( Element::Pointer o1=dynamic_cast<Element*>(&*o) )
124  {
125  el.push_back(FEMP<Element>(o1));
126  continue;
127  }
128  if ( Load::Pointer o1=dynamic_cast<Load*>(&*o) )
129  {
130  load.push_back(FEMP<Load>(o1));
131  continue;
132  }
133 
138  // first we delete the allocated object
139 #ifndef FEM_USE_SMART_POINTERS
140  delete o;
141 #endif
142  o=0;
143 
144  // then we throw an exception
145  throw FEMExceptionIO(__FILE__,__LINE__,"Solver::Read()","Error reading FEM problem stream!");
146 
147  } while ( o );
148 
149 }
150 
151 
155 void Solver::Write( std::ostream& f ) {
156 
157  for(NodeArray::iterator i = node.begin(); i != node.end(); i++)
158  {
159  (*i)->Write(f);
160  }
161  f<<"\n<END> % End of nodes\n\n";
162 
163  for(MaterialArray::iterator i=mat.begin(); i != mat.end(); i++)
164  {
165  (*i)->Write(f);
166  }
167  f<<"\n<END> % End of materials\n\n";
168 
169  for(ElementArray::iterator i=el.begin(); i != el.end(); i++)
170  {
171  (*i)->Write(f);
172  }
173  f<<"\n<END> % End of elements\n\n";
174 
175  for(LoadArray::iterator i=load.begin(); i != load.end(); i++)
176  {
177  (*i)->Write(f);
178  }
179  f<<"\n<END> % End of loads\n\n";
180 
181 }
182 
187 
188  // Clear the list of elements and global freedom numbers in nodes
189  // FIXME: should be removed once Mesh is there
190  for(NodeArray::iterator n=node.begin(); n != node.end(); n++)
191  {
192  (*n)->m_elements.clear();
193  (*n)->ClearDegreesOfFreedom();
194  }
195 
196  for(ElementArray::iterator e=el.begin(); e != el.end(); e++) // step over all elements
197  {
198 
199  // Add the elemens in the nodes list of elements
200  // FIXME: should be removed once Mesh is there
201  unsigned int Npts=(*e)->GetNumberOfNodes();
202  for(unsigned int pt=0; pt<Npts; pt++)
203  {
204  (*e)->GetNode(pt)->m_elements.insert(*e);
205  }
206  }
207 
208 
213  // Start numbering DOFs from 0
214  NGFN=0;
215 
216  // Step over all elements
217  for(ElementArray::iterator e=el.begin(); e != el.end(); e++)
218  {
219  // FIXME: Write a code that checks if two elements are compatible, when they share a node
220  for(unsigned int n=0; n<(*e)->GetNumberOfNodes(); n++)
221  {
222  for(unsigned int dof=0; dof<(*e)->GetNumberOfDegreesOfFreedomPerNode(); dof++)
223  {
224  if( (*e)->GetNode(n)->GetDegreeOfFreedom(dof) == Element::InvalidDegreeOfFreedomID )
225  {
226  (*e)->GetNode(n)->SetDegreeOfFreedom(dof,NGFN);
227  NGFN++;
228  }
229  }
230  }
231  } // end for e
232 
233  // NGFN=Element::GetGlobalDOFCounter()+1;
234  if (NGFN>0) return; // if we got 0 DOF, somebody forgot to define the system...
235 
236 }
237 
242 {
243 
244  // if no DOFs exist in a system, we have nothing to do
245  if (NGFN<=0) return;
246 
247  NMFC=0; // reset number of MFC in a system
248 
254  // search for MFC's in Loads array, because they affect the master stiffness matrix
255  for(LoadArray::iterator l=load.begin(); l != load.end(); l++)
256  {
257  if ( LoadBCMFC::Pointer l1=dynamic_cast<LoadBCMFC*>( &(*(*l))) )
258  {
259  // store the index of an LoadBCMFC object for later
260  l1->Index=NMFC;
261  // increase the number of MFC
262  NMFC++;
263  }
264  }
265 
274 
278  for(ElementArray::iterator e=el.begin(); e != el.end(); e++)
279  {
280  // Call the function that actually moves the element matrix
281  // to the master matrix.
282  this->AssembleElementMatrix(&**e);
283  }
284 
289  for(LoadArray::iterator l2=load.begin(); l2 != load.end(); l2++)
290  {
291  if ( LoadLandmark::Pointer l3=dynamic_cast<LoadLandmark*>( &(*(*l2))) )
292  {
293  l3->AssignToElement(&el);
294  Element::Pointer ep = const_cast<Element*>( l3->el[0] );
295  this->AssembleLandmarkContribution( ep , l3->eta );
296  }
297  }
298 
300 
301 }
302 
304 {
305  // We use LinearSystemWrapper object, to store the K matrix.
306  this->m_ls->SetSystemOrder(N);
307  this->m_ls->InitializeMatrix();
308 }
309 
310 
312 {
313  // Copy the element "landmark" matrix for faster access.
315  e->GetLandmarkContributionMatrix(eta, Le);
316 
317  // ... same for number of DOF
318  int Ne=e->GetNumberOfDegreesOfFreedom();
319 
320  // step over all rows in element matrix
321  for(int j=0; j<Ne; j++)
322  {
323  // step over all columns in element matrix
324  for(int k=0; k<Ne; k++)
325  {
326  // error checking. all GFN should be =>0 and <NGFN
327  if ( e->GetDegreeOfFreedom(j) >= NGFN ||
328  e->GetDegreeOfFreedom(k) >= NGFN )
329  {
330  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::AssembleLandmarkContribution()","Illegal GFN!");
331  }
332 
339  if ( Le[j][k] != Float(0.0) )
340  {
341  this->m_ls->AddMatrixValue( e->GetDegreeOfFreedom(j), e->GetDegreeOfFreedom(k), Le[j][k] );
342  }
343  }
344  }
345 }
346 
348 {
349  // Copy the element stiffness matrix for faster access.
351  e->GetStiffnessMatrix(Ke);
352 
353  // ... same for number of DOF
354  int Ne=e->GetNumberOfDegreesOfFreedom();
355 
356  // step over all rows in element matrix
357  for(int j=0; j<Ne; j++)
358  {
359  // step over all columns in element matrix
360  for(int k=0; k<Ne; k++)
361  {
362  // error checking. all GFN should be =>0 and <NGFN
363  if ( e->GetDegreeOfFreedom(j) >= NGFN ||
364  e->GetDegreeOfFreedom(k) >= NGFN )
365  {
366  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::AssembleElementMatrix()","Illegal GFN!");
367  }
368 
375  if ( Ke[j][k] != Float(0.0) )
376  {
377  this->m_ls->AddMatrixValue( e->GetDegreeOfFreedom(j), e->GetDegreeOfFreedom(k), Ke[j][k] );
378  }
379 
380  }
381 
382  }
383 
384 }
385 
389 void Solver::AssembleF(int dim)
390 {
391 
392  // Vector that stores element nodal loads
394 
395  // Type that stores IDs of fixed DOF together with the values to
396  // which they were fixed.
397  typedef std::map<Element::DegreeOfFreedomIDType,Float> BCTermType;
398  BCTermType bcterm;
399 
400  /* if no DOFs exist in a system, we have nothing to do */
401  if (NGFN<=0) return;
402 
403  /* Initialize the master force vector */
405 
410  for(LoadArray::iterator l=load.begin(); l != load.end(); l++)
411  {
412 
417  Load::Pointer l0=*l;
418 
422  l0->SetSolution(m_ls);
423 
427  if ( LoadNode::Pointer l1=dynamic_cast<LoadNode*>(&*l0) )
428  {
429  // yep, we have a nodal load
430 
431  // size of a force vector in load must match number of DOFs in node
432  if ( (l1->F.size() % l1->m_element->GetNumberOfDegreesOfFreedomPerNode()) != 0 )
433  {
434  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::AssembleF()","Illegal size of a force vector in LoadNode object!");
435  }
436 
437  // we simply copy the load to the force vector
438  for(unsigned int d=0; d < (l1->m_element->GetNumberOfDegreesOfFreedomPerNode()); d++)
439  {
440  Element::DegreeOfFreedomIDType dof=l1->m_element->GetNode(l1->m_pt)->GetDegreeOfFreedom(d);
441  // error checking
442  if ( dof >= NGFN )
443  {
444  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::AssembleF()","Illegal GFN!");
445  }
446 
455  m_ls->AddVectorValue(dof , l1->F[d+l1->m_element->GetNumberOfDegreesOfFreedomPerNode()*dim]);
456  }
457 
458  // that's all there is to DOF loads, go to next load in an array
459  continue;
460  }
461 
462 
466  if ( LoadElement::Pointer l1=dynamic_cast<LoadElement*>(&*l0) )
467  {
468 
469  if ( !(l1->el.empty()) )
470  {
475  for(LoadElement::ElementPointersVectorType::const_iterator i=l1->el.begin(); i != l1->el.end(); i++)
476  {
477 
478  const Element* el0=(*i);
479  // Call the Fe() function of the element that we are applying the load to.
480  // We pass a pointer to the load object as a paramater and a reference to the nodal loads vector.
481  el0->GetLoadVector(Element::LoadPointer(l1),Fe);
482  unsigned int Ne=el0->GetNumberOfDegreesOfFreedom(); // ... element's number of DOF
483  for(unsigned int j=0; j<Ne; j++) // step over all DOF
484  {
485  // error checking
486  if ( el0->GetDegreeOfFreedom(j) >= NGFN )
487  {
488  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::AssembleF()","Illegal GFN!");
489  }
490 
491  // update the master force vector (take care of the correct isotropic dimensions)
492  m_ls->AddVectorValue(el0->GetDegreeOfFreedom(j) , Fe(j+dim*Ne));
493  }
494  }
495 
496  }
497  else
498  {
499 
504  for(ElementArray::iterator e=el.begin(); e != el.end(); e++) // step over all elements in a system
505  {
506  (*e)->GetLoadVector(Element::LoadPointer(l1),Fe); // ... element's force vector
507  unsigned int Ne=(*e)->GetNumberOfDegreesOfFreedom(); // ... element's number of DOF
508 
509  for(unsigned int j=0; j<Ne; j++) // step over all DOF
510  {
511  if ( (*e)->GetDegreeOfFreedom(j) >= NGFN )
512  {
513  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::AssembleF()","Illegal GFN!");
514  }
515 
516  // update the master force vector (take care of the correct isotropic dimensions)
517  m_ls->AddVectorValue((*e)->GetDegreeOfFreedom(j) , Fe(j+dim*Ne));
518 
519  }
520  }
521  }
522 
523  // skip to next load in an array
524  continue;
525  }
526 
530  if ( LoadBCMFC::Pointer l1=dynamic_cast<LoadBCMFC*>(&*l0) )
531  {
532  m_ls->SetVectorValue(NGFN+l1->Index , l1->rhs[dim]);
533 
534  // skip to next load in an array
535  continue;
536  }
537 
541  if ( LoadBC::Pointer l1=dynamic_cast<LoadBC*>(&*l0) )
542  {
543 
544  // Here we just store the values of fixed DOFs. We can't set it here, because
545  // it may be changed by other loads that are applied later.
546  bcterm[ l1->m_element->GetDegreeOfFreedom(l1->m_dof) ]=l1->m_value[dim];
547 
548  // skip to the next load in an array
549  continue;
550  }
551 
558  } // for(LoadArray::iterator l ... )
559 
564  if ( m_ls->IsVectorInitialized(1) )
565  {
566  // Add the vector generated by ApplyBC to the solution vector
567  const unsigned int totGFN=NGFN+NMFC;
568  for( unsigned int i=0; i<totGFN; i++ )
569  {
571  }
572 
573  }
574 
575  // Set the fixed DOFs to proper values
576  for( BCTermType::iterator q=bcterm.begin(); q != bcterm.end(); q++)
577  {
578  m_ls->SetVectorValue(q->first,q->second);
579  }
580 
581 }
582 
587 {
588 }
589 
590 
595 {
596  // Check if master stiffness matrix and master force vector were
597  // properly initialized.
598  if(!m_ls->IsMatrixInitialized())
599  {
600  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::Solve()","Master stiffness matrix was not initialized!");
601  }
602  if(!m_ls->IsVectorInitialized())
603  {
604  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::Solve()","Master force vector was not initialized!");
605  }
606 
607  // Solve the system of linear equations
609  m_ls->Solve();
610 }
611 
617 {
618 }
619 
620 Solver::Float Solver::GetDeformationEnergy(unsigned int SolutionIndex)
621 {
622  Solver::Float U = 0.0f;
623  Element::MatrixType LocalSolution;
624 
625  for(ElementArray::iterator e=el.begin(); e != el.end(); e++)
626  {
627  unsigned int Ne=(*e)->GetNumberOfDegreesOfFreedom();
628  LocalSolution.set_size(Ne,1);
629  // step over all DOFs of element
630  for(unsigned int j=0; j<Ne; j++)
631  {
632  LocalSolution[j][0]=m_ls->GetSolutionValue((*e)->GetDegreeOfFreedom(j),SolutionIndex);
633  }
634 
635  U += (*e)->GetElementDeformationEnergy(LocalSolution);
636  }
637  return U;
638 }
639 
643 void Solver::ApplyBC(int dim, unsigned int matrix)
644 {
645 
646  // Vector with index 1 is used to store force correctios for BCs
647  m_ls->DestroyVector(1);
648 
649  /* Step over all Loads */
650  for(LoadArray::iterator l=load.begin(); l != load.end(); l++)
651  {
652 
657  Load::Pointer l0=*l;
658 
659 
667  if ( LoadBCMFC::Pointer c=dynamic_cast<LoadBCMFC*>(&*l0) )
668  {
669  /* step over all DOFs in MFC */
670  for(LoadBCMFC::LhsType::iterator q=c->lhs.begin();
671  q != c->lhs.end();
672  q++) {
673 
674  /* obtain the GFN of DOF that is in the MFC */
676  q->m_element->GetDegreeOfFreedom(q->dof);
677 
678  /* error checking. all GFN should be =>0 and <NGFN */
679  if ( gfn>=NGFN )
680  {
681  throw FEMExceptionSolution(__FILE__,__LINE__,"Solver::ApplyBC()","Illegal GFN!");
682  }
683 
684  /* set the proper values in matster stiffnes matrix */
685  this->m_ls->SetMatrixValue(gfn, NGFN+c->Index, q->value, matrix);
686  this->m_ls->SetMatrixValue(NGFN+c->Index, gfn, q->value, matrix); // this is a symetric matrix...
687 
688  }
689 
690  // skip to next load in an array
691  continue;
692  }
693 
697  if ( LoadBC::Pointer c=dynamic_cast<LoadBC*>(&*l0) )
698  {
699 
700  Element::DegreeOfFreedomIDType fdof = c->m_element->GetDegreeOfFreedom(c->m_dof);
701  Float fixedvalue=c->m_value[dim];
702 
703 
704  // Copy the corresponding row of the matrix to the vector that will
705  // be later added to the master force vector.
706  // NOTE: We need to copy the whole row first, and then clear it. This
707  // is much more efficient when using sparse matrix storage, than
708  // copying and clearing in one loop.
709 
710  // Get the column indices of the nonzero elements in an array.
712  m_ls->GetColumnsOfNonZeroMatrixElementsInRow(fdof, cols, matrix);
713 
714  // Force vector needs updating only if DOF was not fixed to 0.0.
715  if( fixedvalue != 0.0 )
716  {
717  // Initialize the master force correction vector as required
718  if ( !this->m_ls->IsVectorInitialized(1) )
719  {
720  this->m_ls->InitializeVector(1);
721  }
722 
723  // Step over each nonzero matrix element in a row
724  for(LinearSystemWrapper::ColumnArray::iterator cc=cols.begin(); cc != cols.end(); cc++)
725  {
726  // Get value from the stiffness matrix
727  Float d=this->m_ls->GetMatrixValue(fdof, *cc, matrix);
728 
729  // Store the appropriate value in bc correction vector (-K12*u2)
730  //
731  // See http://titan.colorado.edu/courses.d/IFEM.d/IFEM.Ch04.d/IFEM.Ch04.pdf
732  // chapter 4.1.3 (Matrix Forms of DBC Application Methods) for more info.
733  this->m_ls->AddVectorValue(*cc,-d*fixedvalue,1);
734  }
735  }
736 
737 
738  // Clear that row and column in master matrix
739  for(LinearSystemWrapper::ColumnArray::iterator cc=cols.begin(); cc != cols.end(); cc++)
740  {
741  this->m_ls->SetMatrixValue(fdof,*cc, 0.0, matrix);
742  this->m_ls->SetMatrixValue(*cc,fdof, 0.0, matrix); // this is a symetric matrix
743  }
744  this->m_ls->SetMatrixValue(fdof,fdof, 1.0, matrix); // Set the diagonal element to one
745 
746 
747  // skip to next load in an array
748  continue;
749 
750  }
751 
752  } // end for LoadArray::iterator l
753 }
754 
758 void Solver::InitializeInterpolationGrid(const VectorType& size, const VectorType& bb1, const VectorType& bb2)
759 {
760  // Discard any old image object an create a new one
762 
763  // Set the interpolation grid (image) size, origin and spacing
764  // from the given vectors, so that physical point of v1 is (0,0,0) and
765  // phisical point v2 is (size[0],size[1],size[2]).
766  InterpolationGridType::SizeType image_size={{1,1,1}};
767  for(unsigned int i=0;i<size.size();i++)
768  {
769  image_size[i] = static_cast<InterpolationGridType::SizeType::SizeValueType>( size[i] );
770  }
771  Float image_origin[MaxGridDimensions]={0.0,0.0,0.0};
772  for(unsigned int i=0;i<size.size();i++)
773  {
774  image_origin[i]=bb1[i];
775  }
776  Float image_spacing[MaxGridDimensions]={1.0,1.0,1.0};
777  for(unsigned int i=0;i<size.size();i++)
778  {
779  image_spacing[i]=(bb2[i]-bb1[i])/(image_size[i]-1);
780  }
781 
782  // All regions are the same
783  m_InterpolationGrid->SetRegions(image_size);
784  m_InterpolationGrid->Allocate();
785 
786  // Set origin and spacing
787  m_InterpolationGrid->SetOrigin(image_origin);
788  m_InterpolationGrid->SetSpacing(image_spacing);
789 
790  // Initialize all pointers in interpolation grid image to 0
791  m_InterpolationGrid->FillBuffer(0);
792 
793  VectorType v1,v2;
794 
795  // Fill the interpolation grid with proper pointers to elements
796  for(ElementArray::iterator e=el.begin(); e != el.end(); e++)
797  {
798  // Get square boundary box of an element
799  v1=(*e)->GetNodeCoordinates(0); // lower left corner
800  v2=v1; // upper right corner
801 
802  const unsigned int NumberOfDimensions=(*e)->GetNumberOfSpatialDimensions();
803 
804  for(unsigned int n=1; n < (*e)->GetNumberOfNodes(); n++ )
805  {
806  const VectorType& v=(*e)->GetNodeCoordinates(n);
807  for(unsigned int d=0; d < NumberOfDimensions; d++ )
808  {
809  if( v[d] < v1[d] )
810  {
811  v1[d]=v[d];
812  }
813  if( v[d] > v2[d] )
814  {
815  v2[d]=v[d];
816  }
817  }
818  }
819 
820  // Convert boundary box corner points into discrete image indexes.
822 
824  for(unsigned int i=0;i<MaxGridDimensions;i++)
825  {
826  if ( i < NumberOfDimensions )
827  {
828  vp1[i]=v1[i];
829  vp2[i]=v2[i];
830  }
831  else
832  {
833  vp1[i]=0.0;
834  vp2[i]=0.0;
835  }
836  }
837 
838  // Obtain the Index of BB corner and check whether it is within image.
839  // If it is not, we ignore the entire element.
840  if(!m_InterpolationGrid->TransformPhysicalPointToIndex(vp1,vi1)) continue;
841  if(!m_InterpolationGrid->TransformPhysicalPointToIndex(vp2,vi2)) continue;
842 
844  for( unsigned int i=0; i<MaxGridDimensions; i++ )
845  {
846  region_size[i] = vi2[i]-vi1[i]+1;
847  }
848  InterpolationGridType::RegionType region(vi1,region_size);
849 
850  // Initialize the iterator that will step over all grid points within
851  // element boundary box.
853 
854  //
855  // Update the element pointers in the points defined within the region.
856  //
857  VectorType global_point(NumberOfDimensions); // Point in the image as a vector.
858  VectorType local_point; // Same point in local element coordinate system
859 
860  // Step over all points within the region
861  for(iter.GoToBegin(); !iter.IsAtEnd(); ++iter)
862  {
863  // Note: Iteratior is guarantied to be within image, since the
864  // elements with BB outside are skipped before.
865  m_InterpolationGrid->TransformIndexToPhysicalPoint(iter.GetIndex(),pt);
866  for(unsigned int d=0; d<NumberOfDimensions; d++)
867  {
868  global_point[d]=pt[d];
869  }
870 
871  // If the point is within the element, we update the pointer at
872  // this point in the interpolation grid image.
873  if( (*e)->GetLocalFromGlobalCoordinates(global_point,local_point) )
874  {
875  iter.Set(*e);
876  }
877  } // next point in region
878 
879  } // next element
880 }
881 
882 const Element *
884 {
885  // Add zeros to the end of physical point if necesarry
887  for(unsigned int i=0;i<MaxGridDimensions;i++)
888  {
889  if ( i < pt.size() )
890  {
891  pp[i]=pt[i];
892  }
893  else
894  {
895  pp[i]=0.0;
896  }
897  }
898 
900 
901  // Return value only if given point is within the interpolation grid
902  if( m_InterpolationGrid->TransformPhysicalPointToIndex(pp,index) )
903  {
904  //itk::ContinuousIndex<Float,MaxGridDimensions> ci;
905  //m_InterpolationGrid->TransformPhysicalPointToContinuousIndex(pp,ci);
906  //std::cout<<"In:"<<index<<", ci="<<ci<<", c="<<pp<<"\n";
907  return m_InterpolationGrid->GetPixel(index);
908  }
909  else
910  {
911  // Return 0, if outside the grid.
912  return 0;
913  }
914 }
915 
916 
917 }} // end namespace itk::fem

Generated at Sat May 18 2013 23:38:21 for Orfeo Toolbox with doxygen 1.8.3.1