Orfeo Toolbox  3.16
itkFEMSolverCrankNicolson.cxx
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: $RCSfile: itkFEMSolverCrankNicolson.cxx,v $
5  Language: C++
6  Date: $Date: 2010-03-30 15:20:01 $
7  Version: $Revision: 1.42 $
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 
24 
25 #include "itkFEMLoadNode.h"
26 #include "itkFEMLoadElementBase.h"
27 #include "itkFEMLoadBC.h"
28 #include "itkFEMLoadBCMFC.h"
29 #include "itkFEMLoadLandmark.h"
30 
31 namespace itk {
32 namespace fem {
33 
34 #define TOTE
35 
37 {
52 }
53 
58 {
59 
60  // if no DOFs exist in a system, we have nothing to do
61  if (NGFN<=0) return;
62 
63  Float lhsval;
64  Float rhsval;
65  NMFC=0; // number of MFC in a system
66 
67  // temporary storage for pointers to LoadBCMFC objects
68  typedef std::vector<LoadBCMFC::Pointer> MFCArray;
69  MFCArray mfcLoads;
70 
71  /*
72  * Before we can start the assembly procedure, we need to know,
73  * how many boundary conditions (MFCs) there are in a system.
74  */
75  mfcLoads.clear();
76  // search for MFC's in Loads array, because they affect the master stiffness matrix
77  for(LoadArray::iterator l = load.begin(); l != load.end(); l++)
78  {
79  if ( LoadBCMFC::Pointer l1=dynamic_cast<LoadBCMFC*>( &(*(*l))) )
80  {
81  // store the index of an LoadBCMFC object for later
82  l1->Index=NMFC;
83  // store the pointer to a LoadBCMFC object for later
84  mfcLoads.push_back(l1);
85  // increase the number of MFC
86  NMFC++;
87  }
88  }
89 
95 
96 
100  for(ElementArray::iterator e = el.begin(); e != el.end(); e++)
101  {
102  vnl_matrix<Float> Ke;
103  (*e)->GetStiffnessMatrix(Ke); /*Copy the element stiffness matrix for faster access. */
104  vnl_matrix<Float> Me;
105  (*e)->GetMassMatrix(Me); /*Copy the element mass matrix for faster access. */
106  int Ne=(*e)->GetNumberOfDegreesOfFreedom(); /*... same for element DOF */
107 
108  Me=Me*m_rho;
109 
110  /* step over all rows in in element matrix */
111  for(int j=0; j<Ne; j++)
112  {
113  /* step over all columns in in element matrix */
114  for(int k=0; k<Ne; k++)
115  {
116  /* error checking. all GFN should be =>0 and <NGFN */
117  if ( (*e)->GetDegreeOfFreedom(j) >= NGFN ||
118  (*e)->GetDegreeOfFreedom(k) >= NGFN )
119  {
120  throw FEMExceptionSolution(__FILE__,__LINE__,"SolverCrankNicolson::AssembleKandM()","Illegal GFN!");
121  }
122 
123  /* Here we finaly update the corresponding element
124  * in the master stiffness matrix. We first check if
125  * element in Ke is zero, to prevent zeros from being
126  * allocated in sparse matrix.
127  */
128  if ( Ke(j,k) != Float(0.0) || Me(j,k) != Float(0.0) )
129  {
130  // left hand side matrix
131  lhsval=(Me(j,k) + m_alpha*m_deltaT*Ke(j,k));
132  m_ls->AddMatrixValue( (*e)->GetDegreeOfFreedom(j) ,
133  (*e)->GetDegreeOfFreedom(k),
134  lhsval, SumMatrixIndex );
135  // right hand side matrix
136  rhsval=(Me(j,k) - (1.-m_alpha)*m_deltaT*Ke(j,k));
137  m_ls->AddMatrixValue( (*e)->GetDegreeOfFreedom(j) ,
138  (*e)->GetDegreeOfFreedom(k),
139  rhsval, DifferenceMatrixIndex );
140  }
141  }
142 
143  }
144 
145  }
146 
151  for(LoadArray::iterator l2 = load.begin(); l2 != load.end(); l2++)
152  {
153  if ( LoadLandmark::Pointer l3=dynamic_cast<LoadLandmark*>( &(*(*l2))) )
154  {
155  Element::Pointer ep = const_cast<Element*>( l3->el[0] );
156 
158  ep->GetLandmarkContributionMatrix( l3->eta, Le );
159 
160  int Ne = ep->GetNumberOfDegreesOfFreedom();
161 
162  // step over all rows in element matrix
163  for (int j=0; j<Ne; j++)
164  {
165  // step over all columns in element matrix
166  for (int k=0; k<Ne; k++)
167  {
168  // error checking, all GFN should be >=0 and < NGFN
169  if ( ep->GetDegreeOfFreedom(j) >= NGFN ||
170  ep->GetDegreeOfFreedom(k) >= NGFN )
171  {
172  throw FEMExceptionSolution(__FILE__,__LINE__,"SolverCrankNicolson::AssembleKandM()","Illegal GFN!");
173  }
174 
175  // Now update the corresponding element in the master
176  // stiffness matrix and omit the zeros for the sparseness
177  if ( Le(j,k) != Float(0.0) )
178  {
179  // lhs matrix
180  lhsval = m_alpha*m_deltaT*Le(j,k);
182  ep->GetDegreeOfFreedom(k),
183  lhsval, SumMatrixIndex );
184  //rhs matrix
185  rhsval = (1.-m_alpha)*m_deltaT*Le(j,k);
187  ep->GetDegreeOfFreedom(k),
188  rhsval, DifferenceMatrixIndex );
189  }
190  }
191  }
192  }
193  }
194 
195  /* step over all types of BCs */
196  this->ApplyBC(); // BUG -- are BCs applied appropriately to the problem?
197 }
198 
199 
204 {
205  /* if no DOFs exist in a system, we have nothing to do */
206  if (NGFN<=0) return;
207  AssembleF(dim); // assuming assemblef uses index 0 in vector!
208 
209  typedef std::map<Element::DegreeOfFreedomIDType,Float> BCTermType;
210  BCTermType bcterm;
211 
212  /* Step over all Loads */
213  for(LoadArray::iterator l = load.begin(); l != load.end(); l++)
214  {
215  Load::Pointer l0=*l;
216  if ( LoadBC::Pointer l1=dynamic_cast<LoadBC*>(&*l0) )
217  {
218  bcterm[ l1->m_element->GetDegreeOfFreedom(l1->m_dof) ]=l1->m_value[dim];
219  }
220  } // end for LoadArray::iterator l
221 
222  // Now set the solution t_minus1 vector to fit the BCs
223  for( BCTermType::iterator q = bcterm.begin(); q != bcterm.end(); q++)
224  {
225  m_ls->SetVectorValue(q->first,0.0,SolutionVectorTMinus1Index); //FIXME?
226  m_ls->SetSolutionValue(q->first,0.0,SolutionTMinus1Index); //FIXME?
228  }
229 
232 
233  for (unsigned int index=0; index<NGFN; index++) RecomputeForceVector(index);
234 
235  // Now set the solution and force vector to fit the BCs
236  for( BCTermType::iterator q = bcterm.begin(); q != bcterm.end(); q++)
237  {
238  m_ls->SetVectorValue(q->first,q->second,ForceTIndex);
239  }
240 }
241 
243 {//
244  Float ft = m_ls->GetVectorValue(index,ForceTIndex);
247  Float f=m_deltaT*(m_alpha*ft+(1.-m_alpha)*ftm1)+utm1;
248  m_ls->SetVectorValue(index , f, ForceTIndex);
249 }
250 
255 {
256  /* FIXME - must verify that this is correct use of wrapper */
257  /* FIXME Initialize the solution vector */
259  m_ls->Solve();
260  // call this externally AddToDisplacements();
261 }
262 
263 
265 {
266  // in 1-D domain, we want to find a < b < c , s.t. f(b) < f(a) && f(b) < f(c)
267  // see Numerical Recipes
268 
269  Float Gold=1.618034;
270  Float Glimit=100.0;
271  Float Tiny=1.e-20;
272 
273  Float ax, bx,cx;
274  ax=0.0; bx=1.;
275  Float fc;
276  Float fa=vcl_fabs(EvaluateResidual(ax));
277  Float fb=vcl_fabs(EvaluateResidual(bx));
278 
279  Float ulim,u,r,q,fu,dum;
280 
281  if ( fb > fa )
282  {
283  dum=ax; ax=bx; bx=dum;
284  dum=fb; fb=fa; fa=dum;
285  }
286 
287  cx=bx+Gold*(bx-ax); // first guess for c - the 3rd pt needed to bracket the min
288  fc=vcl_fabs(EvaluateResidual(cx));
289 
290 
291  while (fb > fc /*&& vcl_fabs(ax) < 3. && vcl_fabs(bx) < 3. && vcl_fabs(cx) < 3.*/)
292  {
293  r=(bx-ax)*(fb-fc);
294  q=(bx-cx)*(fb-fa);
295  Float denom=(2.0*GSSign(GSMax(vcl_fabs(q-r),Tiny),q-r));
296  u=(bx)-((bx-cx)*q-(bx-ax)*r)/denom;
297  ulim=bx + Glimit*(cx-bx);
298  if ((bx-u)*(u-cx) > 0.0)
299  {
300  fu=vcl_fabs(EvaluateResidual(u));
301  if (fu < fc)
302  {
303  ax=bx;
304  bx=u;
305  *a=ax; *b=bx; *c=cx;
306  return;
307  }
308  else if (fu > fb)
309  {
310  cx=u;
311  *a=ax; *b=bx; *c=cx;
312  return;
313  }
314 
315  u=cx+Gold*(cx-bx);
316  fu=vcl_fabs(EvaluateResidual(u));
317 
318  }
319  else if ( (cx-u)*(u-ulim) > 0.0)
320  {
321  fu=vcl_fabs(EvaluateResidual(u));
322  if (fu < fc)
323  {
324  bx=cx; cx=u; u=cx+Gold*(cx-bx);
325  fb=fc; fc=fu; fu=vcl_fabs(EvaluateResidual(u));
326  }
327 
328  }
329  else if ( (u-ulim)*(ulim-cx) >= 0.0)
330  {
331  u=ulim;
332  fu=vcl_fabs(EvaluateResidual(u));
333  }
334  else
335  {
336  u=cx+Gold*(cx-bx);
337  fu=vcl_fabs(EvaluateResidual(u));
338  }
339 
340  ax=bx; bx=cx; cx=u;
341  fa=fb; fb=fc; fc=fu;
342 
343  }
344 
345  if ( vcl_fabs(ax) > 1.e3 || vcl_fabs(bx) > 1.e3 || vcl_fabs(cx) > 1.e3)
346  {
347  ax=-2.0; bx=1.0; cx=2.0;
348  } // to avoid crazy numbers caused by bad bracket (u goes nuts)
349 
350  *a=ax; *b=bx; *c=cx;
351 }
352 
354 {
355  // We should now have a, b and c, as well as f(a), f(b), f(c),
356  // where b gives the minimum energy position;
357 
358  Float CGOLD = 0.3819660;
359  Float ZEPS = 1.e-10;
360 
361  Float ax=0.0, bx=1.0, cx=2.0;
362 
363  FindBracketingTriplet(&ax, &bx, &cx);
364 
365  Float xmin;
366 
367  unsigned int iter;
368 
369  Float a,b,d=0.,etemp,fu,fv,fw,fx,p,q,r,tol1,tol2,u,v,w,x,xm;
370 
371  Float e=0.0; // the distance moved on the step before last;
372 
373  a=((ax < cx) ? ax : cx);
374  b=((ax > cx) ? ax : cx);
375 
376  x=w=v=bx;
377  fw=fv=fx=vcl_fabs(EvaluateResidual(x));
378 
379  for (iter = 1; iter <=MaxIters; iter++)
380  {
381  xm=0.5*(a+b);
382  tol2=2.0*(tol1=tol*vcl_fabs(x)+ZEPS);
383  if (vcl_fabs(x-xm) <= (tol2-0.5*(b-a)))
384  {
385  xmin=x;
386  SetEnergyToMin(xmin);
387  return fx;
388  }
389 
390  if (vcl_fabs(e) > tol1)
391  {
392  r=(x-w)*(fx-fv);
393  q=(x-v)*(fx-fw);
394  p=(x-v)*q-(x-w)*r;
395  q=2.0*(q-r);
396  if (q>0.0) p = -1.*p;
397  q=vcl_fabs(q);
398  etemp=e;
399  e=d;
400  if (vcl_fabs(p) >= vcl_fabs(0.5*q*etemp) || p <= q*(a-x) || p >= q*(b-x))
401  d=CGOLD*(e=(x>=xm ? a-x : b-x));
402  else
403  {
404  if (q == 0.0) q=q +ZEPS;
405  d=p/q;
406  u=x+d;
407  if (u-a < tol2 || b-u < tol2) d=GSSign(tol1,xm-x);
408  }
409  }
410  else
411  {
412  d=CGOLD*(e=(x>= xm ? a-x : b-x));
413  }
414 
415  u=(vcl_fabs(d) >= tol1 ? x+d : x + GSSign(tol1,d));
416  fu=vcl_fabs(EvaluateResidual(u));
417  if (fu <= fx)
418  {
419  if ( u >= x ) a=x; else b=x;
420  v=w; w=x;x=u;
421  fv=fw; fw=fx; fx=fu;
422  }
423  else
424  {
425  if (u<x) a = u; else b=u;
426  if (fu <= fw || w ==x)
427  {
428  v=w;
429  w=u;
430  fv=fw;
431  fw=fu;
432  }
433  else if (fu <= fv || v==x || v == w)
434  {
435  v=u;
436  fv=fu;
437  }
438  }
439  }
440  xmin=x;
441  SetEnergyToMin(xmin);
442  return fx;
443 }
444 
446 {
447  // We should now have a, b and c, as well as f(a), f(b), f(c),
448  // where b gives the minimum energy position;
449 
450  Float ax, bx, cx;
451 
452 
453  FindBracketingTriplet(&ax, &bx, &cx);
454  Float xmin,fmin;
455 
456  Float f1,f2,x0,x1,x2,x3;
457 
458  Float R=0.6180339;
459  Float C=(1.0-R);
460 
461  x0=ax;
462  x3=cx;
463  if (vcl_fabs(cx-bx) > vcl_fabs(bx-ax))
464  {
465  x1=bx;
466  x2=bx+C*(cx-bx);
467  }
468  else
469  {
470  x2=bx;
471  x1=bx-C*(bx-ax);
472  }
473  f1=vcl_fabs(EvaluateResidual(x1));
474  f2=vcl_fabs(EvaluateResidual(x2));
475  unsigned int iters=0;
476  while (vcl_fabs(x3-x0) > tol*(vcl_fabs(x1)+vcl_fabs(x2)) && iters < MaxIters)
477  {
478  iters++;
479  if (f2 < f1)
480  {
481  x0=x1; x1=x2; x2=R*x1+C*x3;
482  f1=f2; f2=vcl_fabs(EvaluateResidual(x2));
483  }
484  else
485  {
486  x3=x2; x2=x1; x1=R*x2+C*x0;
487  f2=f1; f1=vcl_fabs(EvaluateResidual(x1));
488  }
489  }
490  if (f1<f2)
491  {
492  xmin=x1;
493  fmin=f1;
494  }
495  else
496  {
497  xmin=x2;
498  fmin=f2;
499  }
500 
501  SetEnergyToMin(xmin);
502  return fmin;
503 }
504 
506 {
507  for (unsigned int j=0; j<NGFN; j++)
508  {
509  Float SolVal;
510  Float FVal;
511 #ifdef LOCE
512  SolVal=xmin*m_ls->GetSolutionValue(j,SolutionTIndex)
514 
515  FVal=xmin*m_ls->GetVectorValue(j,ForceTIndex)
516  +(1.-xmin)*m_ls->GetVectorValue(j,ForceTMinus1Index);
517 #endif
518 #ifdef TOTE
519  SolVal=xmin*m_ls->GetSolutionValue(j,SolutionTIndex);// FOR TOT E
520  FVal=xmin*m_ls->GetVectorValue(j,ForceTIndex);
521 #endif
524  }
525 
526 }
527 
529 {
530  Float DeformationEnergy=0.0;
531  Float iSolVal,jSolVal;
532 
533  for (unsigned int i=0; i<NGFN; i++)
534  {
535 // forming U^T F
536 #ifdef LOCE
537  iSolVal=t*(m_ls->GetSolutionValue(i,SolutionTIndex))
539 #endif
540 #ifdef TOTE
541  iSolVal=t*(m_ls->GetSolutionValue(i,SolutionTIndex))
542  +m_ls->GetSolutionValue(i,TotalSolutionIndex);// FOR TOT E
543 #endif
544 // forming U^T K U
545  Float TempRowVal=0.0;
546  for (unsigned int j=0; j<NGFN; j++)
547  {
548 #ifdef LOCE
549  jSolVal=t*(m_ls->GetSolutionValue(j,SolutionTIndex))
551 #endif
552 #ifdef TOTE
553  jSolVal=t*(m_ls->GetSolutionValue(j,SolutionTIndex))
554  +m_ls->GetSolutionValue(j,TotalSolutionIndex);// FOR TOT E
555 #endif
556  TempRowVal += m_ls->GetMatrixValue(i,j,SumMatrixIndex)*jSolVal;
557  }
558  DeformationEnergy += iSolVal*TempRowVal;
559  }
560  return DeformationEnergy;
561 }
562 
563 
565 {
566 
567  Float ForceEnergy=0.0,FVal=0.0;
568  Float DeformationEnergy=0.0;
569  Float iSolVal,jSolVal;
570 
571  for (unsigned int i=0; i<NGFN; i++)
572  {
573 // forming U^T F
574 #ifdef LOCE
575  iSolVal=t*(m_ls->GetSolutionValue(i,SolutionTIndex))
578  FVal=t*FVal+(1.-t)*m_ls->GetVectorValue(i,ForceTMinus1Index);
579 
580  ForceEnergy += iSolVal*FVal;
581 #endif
582 #ifdef TOTE
583  FVal=FVal+0.0;
584  iSolVal=t*(m_ls->GetSolutionValue(i,SolutionTIndex))
585  +m_ls->GetSolutionValue(i,TotalSolutionIndex);// FOR TOT E
586  ForceEnergy += iSolVal*(m_ls->GetVectorValue(i,ForceTotalIndex)+
587  t*m_ls->GetVectorValue(i,ForceTIndex));// FOR TOT E
588 #endif
589 // forming U^T K U
590  Float TempRowVal=0.0;
591  for (unsigned int j=0; j<NGFN; j++)
592  {
593 #ifdef LOCE
594  jSolVal=t*(m_ls->GetSolutionValue(j,SolutionTIndex))
596 #endif
597 #ifdef TOTE
598  jSolVal=t*(m_ls->GetSolutionValue(j,SolutionTIndex))
599  +m_ls->GetSolutionValue(j,TotalSolutionIndex);// FOR TOT E
600 #endif
601  TempRowVal += m_ls->GetMatrixValue(i,j,SumMatrixIndex)*jSolVal;
602  }
603  DeformationEnergy += iSolVal*TempRowVal;
604  }
605  Float Energy=(Float) vcl_fabs(DeformationEnergy-ForceEnergy);
606  return Energy;
607 }
608 
614 {
619  Float maxs=0.0,CurrentTotSolution,CurrentSolution,CurrentForce;
620  Float mins2=0.0, maxs2=0.0;
621  Float absmax=0.0;
622 
623  for(unsigned int i=0;i<NGFN;i++)
624  {
625 #ifdef TOTE
626  CurrentSolution=m_ls->GetSolutionValue(i,SolutionTIndex);
627 #endif
628  if (CurrentSolution < mins2 )
629  {
630  mins2=CurrentSolution;
631  }
632  else if (CurrentSolution > maxs2 )
633  {
634  maxs2=CurrentSolution;
635  }
636  if (vcl_fabs(CurrentSolution) > absmax) absmax=vcl_fabs(CurrentSolution);
637 
638 // note: set rather than add - i.e. last solution of system not total solution
639 #ifdef LOCE
640  CurrentSolution=optimum*m_ls->GetSolutionValue(i,SolutionTIndex)
642  CurrentForce=optimum*m_ls->GetVectorValue(i,ForceTIndex)
643  +(1.-optimum)*m_ls->GetVectorValue(i,ForceTMinus1Index);
644  m_ls->SetVectorValue(i,CurrentSolution,SolutionVectorTMinus1Index);
645  m_ls->SetSolutionValue(i,CurrentSolution,SolutionTMinus1Index);
646  m_ls->SetVectorValue(i , CurrentForce, ForceTMinus1Index); // now set t minus one force vector correctly
647 #endif
648 #ifdef TOTE
649  CurrentSolution=optimum*CurrentSolution;
650  CurrentForce=optimum*m_ls->GetVectorValue(i,ForceTIndex);
651  m_ls->SetVectorValue(i,CurrentSolution,SolutionVectorTMinus1Index); // FOR TOT E
652  m_ls->SetSolutionValue(i,CurrentSolution,SolutionTMinus1Index); // FOR TOT E
653  m_ls->SetVectorValue(i,CurrentForce,ForceTMinus1Index);
654 #endif
655 
656  m_ls->AddSolutionValue(i,CurrentSolution,TotalSolutionIndex);
657  m_ls->AddVectorValue(i , CurrentForce, ForceTotalIndex);
658  CurrentTotSolution=m_ls->GetSolutionValue(i,TotalSolutionIndex);
659 
660  if ( vcl_fabs(CurrentTotSolution) > maxs )
661  {
662  maxs=vcl_fabs(CurrentTotSolution);
663  }
664 
665 
666  }
667 
668  m_CurrentMaxSolution=absmax;
669 }
670 
675 {
680  Float mins=0.0, maxs=0.0;
681  Float mins2=0.0, maxs2=0.0;
682  for(unsigned int i=0;i<NGFN;i++)
683  {
684  Float CurrentSolution=m_ls->GetSolutionValue(i,SolutionTIndex);
685  if (CurrentSolution < mins2 ) mins2=CurrentSolution;
686  else if (CurrentSolution > maxs2 ) maxs2=CurrentSolution;
687  CurrentSolution=m_ls->GetSolutionValue(i,TotalSolutionIndex);
688  if (CurrentSolution < mins ) mins=CurrentSolution;
689  else if (CurrentSolution > maxs ) maxs=CurrentSolution;
690  }
691 }
692 
698 {
699 
700  Float maxs=0.0;
701  for(unsigned int i=0;i<NGFN;i++)
702  {
705  Float newsol=t*(temp)+(1.-t)*temp2;
709  if ( newsol > maxs ) maxs=newsol;
710  }
711 }
712 
714 {
715  for(unsigned int i=0;i<NGFN;i++)
716  {
717  m_ls->SetVectorValue(i,0.0,which);
718  }
719 }
720 
722 {
723  std::cout << " printing current displacements " << std::endl;
724  for(unsigned int i=0;i<NGFN;i++)
725  {
726  std::cout << m_ls->GetSolutionValue(i,TotalSolutionIndex) << std::endl;
727  }
728 }
729 
731 {
732  std::cout << " printing current forces " << std::endl;
733  for(unsigned int i=0;i<NGFN;i++)
734  {
735  std::cout << m_ls->GetVectorValue(i,ForceTIndex) << std::endl;
736  }
737 }
738 
739 }} // end namespace itk::fem

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