Orfeo Toolbox  3.16
itkFEMElementBase.cxx
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: $RCSfile: itkFEMElementBase.cxx,v $
5  Language: C++
6  Date: $Date: 2009-01-29 20:09:13 $
7  Version: $Revision: 1.39 $
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 "itkFEMElementBase.h"
24 #include "vnl/algo/vnl_svd.h"
25 #include "vnl/algo/vnl_qr.h"
26 
27 namespace itk {
28 namespace fem {
29 
30 #ifdef FEM_BUILD_VISUALIZATION
31 
33 double Element::DC_Scale=1000.0;
34 
36 double &Element::Node::DC_Scale=Element::DC_Scale;
37 
38 
42 void Element::Node::Draw(CDC* pDC, Solution::ConstPointer sol) const
43 {
44  // We can only draw 2D nodes here
45  if(m_coordinates.size() != 2)
46  {
47  return;
48  }
49 
50  // Normally we draw a white circle.
51  CPen pen(PS_SOLID, 0, (COLORREF) RGB(0,0,0) );
52  CBrush brush( RGB(255,255,255) );
53 
54  CPen* pOldpen=pDC->SelectObject(&pen);
55  CBrush* pOldbrush=pDC->SelectObject(&brush);
56 
57  int x1=m_coordinates[0]*DC_Scale;
58  int y1=m_coordinates[1]*DC_Scale;
59  x1 += sol->GetSolutionValue(this->GetDegreeOfFreedom(0))*DC_Scale;
60  y1 += sol->GetSolutionValue(this->GetDegreeOfFreedom(1))*DC_Scale;
61 
62  CPoint r1=CPoint(0,0);
63  CPoint r=CPoint(5,5);
64 
65  pDC->DPtoLP(&r1);
66  pDC->DPtoLP(&r);
67  r=r-r1;
68 
69  pDC->Ellipse(x1-r.x, y1-r.y, x1+r.x, y1+r.y);
70 
71  pDC->SelectObject(pOldbrush);
72  pDC->SelectObject(pOldpen);
73 
74 }
75 
76 #endif
77 
81 void Element::Node::Read( std::istream& f, void* info )
82 {
83  unsigned int n;
84 
88  Superclass::Read(f,info);
89 
90  /*
91  * Read and set node coordinates
92  */
93  this->SkipWhiteSpace(f); f>>n; if(!f) goto out;
94  this->m_coordinates.set_size(n);
95  this->SkipWhiteSpace(f); f>>this->m_coordinates; if(!f) goto out;
96 
97 out:
98 
99  if( !f )
100  {
101  throw FEMExceptionIO(__FILE__,__LINE__,"Element::Node::Read()","Error reading FEM node!");
102  }
103 }
104 
105 /*
106  * Write the Node to the output stream
107  */
108 void Element::Node::Write( std::ostream& f ) const
109 {
114 
119  /* write the value of dof */
120  f<<"\t"<<this->m_coordinates.size();
121  f<<" "<<this->m_coordinates<<"\t% Node coordinates"<<"\n";
122 
124  if (!f)
125  {
126  throw FEMExceptionIO(__FILE__,__LINE__,"Element::Node::Write()","Error writing FEM node!");
127  }
128 }
129 
131 
135 /*
136  * Compute and return element stiffnes matrix (Ke) in global coordinate
137  * system.
138  * The base class provides a general implementation which only computes
139  *
140  * b T
141  * int B(x) D B(x) dx
142  * a
143  *
144  * using the Gaussian numeric integration method. The function calls
145  * GetIntegrationPointAndWeight() / GetNumberOfIntegrationPoints() to obtain
146  * the integration points. It also calls the GetStrainDisplacementMatrix()
147  * and GetMaterialMatrix() member functions.
148  *
149  * \param Ke Reference to the resulting stiffnes matrix.
150  *
151  * \note This is a very generic implementation of the stiffness matrix
152  * that is suitable for any problem/element definition. A specifc
153  * element may override this implementation with its own simple one.
154  */
156 {
157  // B and D matrices
158  MatrixType B,D;
159  this->GetMaterialMatrix( D );
160 
161  unsigned int Nip=this->GetNumberOfIntegrationPoints();
162 
163  VectorType ip;
164  Float w;
165  MatrixType J;
166  MatrixType shapeDgl;
167  MatrixType shapeD;
168 
169  // Calculate the contribution of 1st int. point to initialize
170  // the Ke matrix to proper number of elements.
171  this->GetIntegrationPointAndWeight(0,ip,w);
172  this->ShapeFunctionDerivatives(ip,shapeD);
173  this->Jacobian(ip,J,&shapeD);
174  this->ShapeFunctionGlobalDerivatives(ip,shapeDgl,&J,&shapeD);
175 
176  this->GetStrainDisplacementMatrix( B, shapeDgl );
177  Float detJ=this->JacobianDeterminant( ip, &J );
178  Ke=detJ*w*B.transpose()*D*B; // FIXME: write a more efficient way of computing this.
179 
180  // Add contributions of other int. points to the Ke
181  for(unsigned int i=1; i<Nip; i++)
182  {
183  this->GetIntegrationPointAndWeight(i,ip,w);
184  this->ShapeFunctionDerivatives(ip,shapeD);
185  this->Jacobian(ip,J,&shapeD);
186  this->ShapeFunctionGlobalDerivatives(ip,shapeDgl,&J,&shapeD);
187 
188  this->GetStrainDisplacementMatrix( B, shapeDgl );
189  detJ=this->JacobianDeterminant( ip, &J );
190  Ke += detJ*w*B.transpose()*D*B; // FIXME: write a more efficient way of computing this.
191  }
192 }
193 
194 Element::VectorType Element::GetStrainsAtPoint(const VectorType& pt, const Solution& sol, unsigned int index) const
195 // NOTE: pt should be in local coordinates already
196 {
197  MatrixType B;
198  VectorType e, u;
199  MatrixType J, shapeD, shapeDgl;
200 
201  this->ShapeFunctionDerivatives(pt, shapeD);
202  this->Jacobian(pt, J, &shapeD);
203  this->ShapeFunctionGlobalDerivatives(pt, shapeDgl, &J, &shapeD);
204  this->GetStrainDisplacementMatrix(B, shapeDgl);
205 
206  u = this->InterpolateSolution(pt, sol, index);
207 
208  e = B*u;
209 
210  return e;
211 }
212 
214  const VectorType& e,
215  const Solution& itkNotUsed(sol),
216  unsigned int itkNotUsed(index)) const
217 // NOTE: pt should be in local coordinates already
218 {
219  MatrixType D;
220  VectorType sigma;
221 
222  this->GetMaterialMatrix(D);
223 
224  sigma = D*e;
225 
226  return sigma;
227 }
228 
230 {
231  // Provides the contribution of a landmark to the element stiffness
232  // matrix
234 
235  const unsigned int NnDOF=this->GetNumberOfDegreesOfFreedomPerNode();
236  const unsigned int Nnodes=this->GetNumberOfNodes();
237  const unsigned int NDOF = GetNumberOfDegreesOfFreedom();
238  const unsigned int Nip=this->GetNumberOfIntegrationPoints(0);
239 
240  Le.set_size(NDOF,NDOF); // resize the target matrix object
241  Le.fill(0.0);
242 
243  Float w;
244  VectorType ip, shape;
245 
246  for(unsigned int i=0; i<Nip; i++)
247  {
248  this->GetIntegrationPointAndWeight(i,ip,w,0);
249  shape=this->ShapeFunctions(ip);
250 
251  for(unsigned int ni=0; ni<Nnodes; ni++)
252  {
253  for(unsigned int nj=0; nj<Nnodes; nj++)
254  {
255  Float m=w*shape[ni]*shape[nj];
256  for(unsigned int d=0; d<NnDOF; d++)
257  {
258  Le[ni*NnDOF+d][nj*NnDOF+d] += m;
259  }
260  }
261  }
262  }
263 
264  Le = Le / (eta );
265 }
266 
267 
269 {
270 
271  MatrixType U;
272 
274  this->GetStiffnessMatrix(Ke);
275 
276  U=LocalSolution.transpose()*Ke*LocalSolution;
277 
278  return U[0][0];
279 }
280 
282 {
283  /*
284  * If the function is not overiden, we compute consistent mass matrix
285  * by integrating the shape functions over the element domain. The element
286  * density is assumed one. If this is not the case, the GetMassMatrix in a
287  * derived class must be overriden and the Me matrix corrected accordingly.
288  */
290 
291  const unsigned int NnDOF=this->GetNumberOfDegreesOfFreedomPerNode();
292  const unsigned int Nnodes=this->GetNumberOfNodes();
293  const unsigned int NDOF = GetNumberOfDegreesOfFreedom();
294  const unsigned int Nip=this->GetNumberOfIntegrationPoints(0);
295 
296  Me.set_size(NDOF,NDOF); // resize the target matrix object
297  Me.fill(0.0);
298 
299  Float w;
300  VectorType ip, shape;
301  MatrixType J, shapeD;
302 
303  for(unsigned int i=0; i<Nip; i++)
304  {
305  this->GetIntegrationPointAndWeight(i,ip,w,0);
306  shape=this->ShapeFunctions(ip);
307  this->ShapeFunctionDerivatives(ip,shapeD);
308  this->Jacobian(ip,J,&shapeD);
309  Float detJ=this->JacobianDeterminant( ip, &J );
310 
311  for(unsigned int ni=0; ni<Nnodes; ni++)
312  {
313  for(unsigned int nj=0; nj<Nnodes; nj++)
314  {
315  Float m=detJ*w*shape[ni]*shape[nj];
316  for(unsigned int d=0; d<NnDOF; d++)
317  {
318  Me[ni*NnDOF+d][nj*NnDOF+d] += m;
319  }
320  }
321  }
322  }
323 
324 }
325 
327 Element::InterpolateSolution( const VectorType& pt, const Solution& sol, unsigned int solutionIndex ) const
328 {
329 
331  VectorType shapef = this->ShapeFunctions(pt);
332  Float value;
333 
334  const unsigned int Nnodes=this->GetNumberOfNodes();
335  const unsigned int Ndofs_per_node=this->GetNumberOfDegreesOfFreedomPerNode();
336 
337  for(unsigned int f=0; f<Ndofs_per_node; f++)
338  {
339  value=0.0;
340 
341  for(unsigned int n=0; n<Nnodes; n++)
342  {
343  value += shapef[n] * sol.GetSolutionValue( this->GetNode(n)->GetDegreeOfFreedom(f) , solutionIndex);
344  }
345 
346  vec[f]=value;
347 
348  }
349 
350  return vec;
351 
352 }
353 
355 Element::InterpolateSolutionN( const VectorType& pt, const Solution& sol, unsigned int f , unsigned int solutionIndex ) const
356 {
357 
358  Float value=0.0;
359 
360  VectorType shapef = this->ShapeFunctions(pt);
361  unsigned int Nnodes=this->GetNumberOfNodes();
362  for(unsigned int n=0; n<Nnodes; n++)
363  {
364  value += shapef[n] * sol.GetSolutionValue( this->GetNode(n)->GetDegreeOfFreedom(f), solutionIndex );
365  }
366  return value;
367 
368 }
369 
371 
375 void
376 Element::Jacobian( const VectorType& pt, MatrixType& J, const MatrixType* pshapeD ) const
377 {
378  MatrixType* pshapeDlocal=0;
379 
380  // If derivatives of shape functions were not provided, we
381  // need to compute them here
382  if(pshapeD==0)
383  {
384  pshapeDlocal=new MatrixType();
385  this->ShapeFunctionDerivatives( pt, *pshapeDlocal );
386  pshapeD=pshapeDlocal;
387  }
388 
389  const unsigned int Nn=pshapeD->columns();
390  const unsigned int Ndims=this->GetNumberOfSpatialDimensions();
391 
392  MatrixType coords(Nn, Ndims);
393 
394  for( unsigned int n=0; n<Nn; n++ )
395  {
396  VectorType p=this->GetNodeCoordinates(n);
397  coords.set_row(n,p);
398  }
399 
400  J=(*pshapeD)*coords;
401 
402  // Destroy local copy of derivatives of shape functions, if
403  // they were computed.
404  delete pshapeDlocal;
405 
406 }
407 
410 {
411  MatrixType* pJlocal=0;
412 
413  // If Jacobian was not provided, we
414  // need to compute it here
415  if(pJ==0)
416  {
417  pJlocal=new MatrixType();
418  this->Jacobian( pt, *pJlocal );
419  pJ=pJlocal;
420  }
421 
422  // Float det=vnl_svd<Float>(*pJ).determinant_magnitude();
423  Float det=vnl_qr<Float>(*pJ).determinant();
424 
425  delete pJlocal;
426 
427  return det;
428 
429 }
430 
431 void
432 Element::JacobianInverse( const VectorType& pt, MatrixType& invJ, const MatrixType* pJ ) const
433 {
434 
435  MatrixType* pJlocal=0;
436 
437  // If Jacobian was not provided, we
438  // need to compute it here
439  if(pJ==0)
440  {
441  pJlocal=new MatrixType();
442  this->Jacobian( pt, *pJlocal );
443  pJ=pJlocal;
444  }
445 
446  // invJ=vnl_svd_inverse<Float>(*pJ);
447  invJ=vnl_qr<Float>(*pJ).inverse();
448 
449  delete pJlocal;
450 
451 }
452 
453 void Element::ShapeFunctionGlobalDerivatives( const VectorType& pt, MatrixType& shapeDgl, const MatrixType* pJ, const MatrixType* pshapeD ) const
454 {
455 
456  MatrixType* pshapeDlocal=0;
457  MatrixType* pJlocal=0;
458 
459  // If derivatives of shape functions were not provided, we
460  // need to compute them here
461  if(pshapeD==0)
462  {
463  pshapeDlocal=new MatrixType();
464  this->ShapeFunctionDerivatives( pt, *pshapeDlocal );
465  pshapeD=pshapeDlocal;
466  }
467 
468  // If Jacobian was not provided, we
469  // need to compute it here
470  if(pJ==0)
471  {
472  pJlocal=new MatrixType();
473  this->Jacobian( pt, *pJlocal, pshapeD );
474  pJ=pJlocal;
475  }
476 
477  MatrixType invJ;
478  this->JacobianInverse( pt, invJ, pJ );
479 
480  shapeDgl=invJ*(*pshapeD);
481 
482  delete pJlocal;
483  delete pshapeDlocal;
484 
485 }
486 
489 {
490  unsigned int Nnodes=this->GetNumberOfNodes();
491  MatrixType nc(this->GetNumberOfSpatialDimensions(),Nnodes);
492  for(unsigned int n=0; n<Nnodes; n++)
493  {
494  nc.set_column( n,this->GetNodeCoordinates(n) );
495  }
496 
497  VectorType shapeF = ShapeFunctions(pt);
498 
499  return nc*shapeF;
500 
501 }
502 
503 // Gauss-Legendre integration rule constants
505 {
506  { 0.0 },
507  { 0.000000000000000 },
508  { 0.577350269189626,-0.577350269189626 },
509  { 0.774596669241483, 0.000000000000000,-0.774596669241483 },
510  { 0.861136311594053, 0.339981043584856,-0.339981043584856,-0.861136311594053 },
511  { 0.906179845938664, 0.538469310105683, 0.000000000000000,-0.538469310105683,-0.906179845938664},
512  { 0.932469514203152, 0.661209386466264, 0.238619186083197,-0.238619186083197,-0.661209386466264,-0.932469514203152 },
513  { 0.949107912342759, 0.741531185599394, 0.405845151377397, 0.000000000000000,-0.405845151377397,-0.741531185599394,-0.949107912342759 },
514  { 0.960289856497536, 0.796666477413627, 0.525532409916329, 0.183434642495650,-0.183434642495650,-0.525532409916329,-0.796666477413627,-0.960289856497536 },
515  { 0.968160239507626, 0.836031107326636, 0.613371432700590, 0.324253423403809, 0.000000000000000,-0.324253423403809,-0.613371432700590,-0.836031107326636,-0.968160239507626 },
516  { 0.973906528517172, 0.865063366688985, 0.679409568299024, 0.433395394129247, 0.148874338981631,-0.148874338981631,-0.433395394129247,-0.679409568299024,-0.865063366688985,-0.973906528517172 }
517 };
518 
520 {
521  { 0.0 },
522  { 2.000000000000000 },
523  { 1.000000000000000, 1.000000000000000 },
524  { 0.555555555555555, 0.888888888888889, 0.555555555555555 },
525  { 0.347854845137454, 0.652145154862546, 0.652145154862546, 0.347854845137454 },
526  { 0.236926885056189, 0.478628670499366, 0.568888888888889, 0.478628670499366, 0.236926885056189 },
527  { 0.171324492379170, 0.360761573048139, 0.467913934572691, 0.467913934572691, 0.360761573048139, 0.171324492379170 },
528  { 0.129484966168869, 0.279705391489277, 0.381830050505119, 0.417959183673469, 0.381830050505119, 0.279705391489277, 0.129484966168869 },
529  { 0.101228536290376, 0.222381034453374, 0.313706645877887, 0.362683783378362, 0.362683783378362, 0.313706645877887, 0.222381034453374, 0.101228536290376 },
530  { 0.081274388361575, 0.180648160694858, 0.260610696402935, 0.312347077040003, 0.330239355001260, 0.312347077040003, 0.260610696402935, 0.180648160694858, 0.081274388361575 },
531  { 0.066671344308688, 0.149451349150581, 0.219086362515982, 0.269266719309996, 0.295524224714753, 0.295524224714753, 0.269266719309996, 0.219086362515982, 0.149451349150581, 0.066671344308688 }
532 };
533 
534 // Register Node class with FEMObjectFactory
536 
537 }} // end namespace itk::fem

Generated at Sat Feb 2 2013 23:36:42 for Orfeo Toolbox with doxygen 1.8.1.1