LifeV
lifev/structure/examples/example_bodyForces/main.cpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 /**
27  \file main.cpp
28 
29  This test is the case of traction of a cube. It does not use the symmetry BCs
30  This test uses the FESpace which is the standard in LifeV and the ETFESpace
31  The FESpace is used for the BCs of Neumann type since in the ET branch there
32  is not the integration over the boundary faces.
33 
34  \author Paolo Tricerri <paolo.tricerri@epfl.ch>
35  \date 1861-03-17
36  */
37 
38 #ifdef TWODIM
39 #error test_structure cannot be compiled in 2D
40 #endif
41 
42 // Tell the compiler to ignore specific kind of warnings:
43 #pragma GCC diagnostic ignored "-Wunused-variable"
44 #pragma GCC diagnostic ignored "-Wunused-parameter"
45 
46 #include <Epetra_ConfigDefs.h>
47 #ifdef EPETRA_MPI
48 #include <mpi.h>
49 #include <Epetra_MpiComm.h>
50 #else
51 #include <Epetra_SerialComm.h>
52 #endif
53 
54 //Tell the compiler to restore the warning previously silented
55 #pragma GCC diagnostic warning "-Wunused-variable"
56 #pragma GCC diagnostic warning "-Wunused-parameter"
57 
58 #include <lifev/core/LifeV.hpp>
59 #include <lifev/core/algorithm/PreconditionerIfpack.hpp>
60 #include <lifev/core/algorithm/PreconditionerML.hpp>
61 
62 #include <lifev/core/array/MapEpetra.hpp>
63 
64 #include <lifev/core/fem/TimeAdvance.hpp>
65 #include <lifev/core/fem/TimeAdvanceBDF.hpp>
66 
67 #include <lifev/core/mesh/MeshData.hpp>
68 #include <lifev/core/mesh/MeshPartitioner.hpp>
69 
70 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp>
71 
72 #include <lifev/structure/solver/StructuralConstitutiveLaw.hpp>
73 #include <lifev/structure/solver/StructuralOperator.hpp>
74 #include <lifev/structure/solver/isotropic/NeoHookeanMaterialNonLinear.hpp>
75 #include <lifev/structure/solver/anisotropic/StructuralAnisotropicConstitutiveLaw.hpp>
76 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp>
77 
78 
79 #include <lifev/core/filter/ExporterEnsight.hpp>
80 #ifdef HAVE_HDF5
81 #include <lifev/core/filter/ExporterHDF5.hpp>
82 #endif
83 #include <lifev/core/filter/ExporterEmpty.hpp>
84 
85 //Includes for the Expression Template
86 #include <lifev/eta/fem/ETFESpace.hpp>
87 
88 #include <iostream>
89 
90 #include "ud_functions.hpp"
91 
92 
93 using namespace LifeV;
94 
95 int returnValue = EXIT_FAILURE;
97 
98 namespace
99 {
100 static bool regIF = (PRECFactory::instance().registerProduct ( "Ifpack", &createIfpack ) );
101 static bool regML = (PRECFactory::instance().registerProduct ( "ML", &createML ) );
102 }
103 
104 
106 {
107  std::string stringList = list;
108  std::set<UInt> setList;
109  if ( list == "" )
110  {
111  return setList;
112  }
113  size_t commaPos = 0;
114  while ( commaPos != std::string::npos )
115  {
116  commaPos = stringList.find ( "," );
117  setList.insert ( atoi ( stringList.substr ( 0, commaPos ).c_str() ) );
118  stringList = stringList.substr ( commaPos + 1 );
119  }
120  setList.insert ( atoi ( stringList.c_str() ) );
121  return setList;
122 }
123 
124 
125 class Structure
126 {
127 public:
128 
129  // Public typedefs
135 
138 
139  // typedefs for fibers
140  // std function for fiber direction
141  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > analyticalFunction_Type;
143 
144  typedef std::function<VectorSmall<3> ( Real const&, Real const&, Real const&, Real const& ) > bodyFunction_Type;
146 
147 
150 
152 
153  /** @name Constructors, destructor
154  */
155  //@{
156  Structure ( int argc,
157  char** argv,
158  std::shared_ptr<Epetra_Comm> structComm );
159 
161  {}
162  //@}
163 
164  //@{
165  void run()
166  {
167  run3d();
168  }
169  //@}
170 
171 protected:
172 
173 private:
174 
175  /**
176  * run the 2D driven cylinder simulation
177  */
178  void run2d();
179 
180  /**
181  * run the 3D driven cylinder simulation
182  */
183  void run3d();
184 
185 private:
186  struct Private;
188 };
189 
190 
191 
192 struct Structure::Private
193 {
195  rho (1), poisson (1), young (1), bulk (1), alpha (1), gamma (1)
196  {}
197  double rho, poisson, young, bulk, alpha, gamma;
198 
200 
202 };
203 
204 
205 
206 Structure::Structure ( int argc,
207  char** argv,
208  std::shared_ptr<Epetra_Comm> structComm) :
209  parameters ( new Private() )
210 {
211  GetPot command_line (argc, argv);
212  std::string data_file_name = command_line.follow ("data", 2, "-f", "--file");
213  GetPot dataFile ( data_file_name );
214  parameters->data_file_name = data_file_name;
215 
216  parameters->comm = structComm;
217  int ntasks = parameters->comm->NumProc();
218 
219  if (!parameters->comm->MyPID() )
220  {
221  std::cout << "My PID = " << parameters->comm->MyPID() << " out of " << ntasks << " running." << std::endl;
222  }
223 }
224 
225 
226 
227 void
228 Structure::run2d()
229 {
230  std::cout << "2D cylinder test case is not available yet\n";
231 }
232 
233 
234 
235 void
236 Structure::run3d()
237 {
238  bool verbose = (parameters->comm->MyPID() == 0);
239 
240  //! Number of boundary conditions for the velocity and mesh motion
241  std::shared_ptr<BCHandler> BCh ( new BCHandler() );
242 
243  //! dataElasticStructure
244  GetPot dataFile ( parameters->data_file_name.c_str() );
245 
246  std::shared_ptr<StructuralConstitutiveLawData> dataStructure (new StructuralConstitutiveLawData( ) );
247  dataStructure->setup (dataFile);
248 
249  dataStructure->showMe();
250 
251  MeshData meshData;
252  meshData.setup (dataFile, "solid/space_discretization");
253 
254  std::shared_ptr<RegionMesh<LinearTetra> > fullMeshPtr (new RegionMesh<LinearTetra> ( parameters->comm ) );
255  readMesh (*fullMeshPtr, meshData);
256 
257  MeshPartitioner< RegionMesh<LinearTetra> > meshPart ( fullMeshPtr, parameters->comm );
258  std::shared_ptr<RegionMesh<LinearTetra> > localMeshPtr (new RegionMesh<LinearTetra> ( parameters->comm ) );
259  localMeshPtr = meshPart.meshPartition();
260 
261  std::string dOrder = dataFile ( "solid/space_discretization/order", "P1");
262 
263  //Mainly used for BCs assembling (Neumann type)
264  solidFESpacePtr_Type dFESpace ( new solidFESpace_Type (localMeshPtr, dOrder, 3, parameters->comm) );
265  solidETFESpacePtr_Type dETFESpace ( new solidETFESpace_Type (meshPart, & (dFESpace->refFE() ), & (dFESpace->fe().geoMap() ), parameters->comm) );
266 
267  vectorFiberFunctionPtr_Type pointerToVectorOfFamilies;
268  listOfFiberDirections_Type fiberDirections;
269  fibersDirectionList setOfFiberFunctions;
270 
271  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
272  {
273  // Setting the fibers
274  pointerToVectorOfFamilies.reset( new vectorFiberFunction_Type( ) );
275  (*pointerToVectorOfFamilies).resize( dataStructure->numberFibersFamilies( ) );
276 
277  fiberDirections.resize( dataStructure->numberFibersFamilies( ) );
278 
279  // if ( verbose )
280  // {
281  // std::cout << "Size of the number of families: " << (*pointerToVectorOfFamilies).size() << std::endl;
282  // }
283 
284  setOfFiberFunctions.setupFiberDefinitions( dataStructure->numberFibersFamilies( ) );
285  }
286 
287  std::string timeAdvanceMethod = dataFile ( "solid/time_discretization/method", "BDF");
288 
289  timeAdvance_Type timeAdvance ( TimeAdvanceFactory::instance().createObject ( timeAdvanceMethod ) );
290 
291  UInt OrderDev = 2;
292 
293  timeAdvance->setup (dataStructure->dataTimeAdvance()->orderBDF() , OrderDev);
294 
295  timeAdvance->setTimeStep (dataStructure->dataTime()->timeStep() );
296  //timeAdvance->showMe();
297 
298  //! #################################################################################
299  //! BOUNDARY CONDITIONS + SOURCE TERMS
300  //! #################################################################################
301  std::vector <ID> compx (1), compy (1), compz (1), compxy (2), compxz (2), compyz (2);
302  compx[0] = 0;
303  compy[0] = 1, compz[0] = 2;
304  compxy[0] = 0;
305  compxy[1] = 1;
306  compxz[0] = 0;
307  compxz[1] = 2;
308  compyz[0] = 1;
309  compyz[1] = 2;
310 
311  BCFunctionBase zero (bcZero);
312  BCFunctionBase nonZero;
313  BCFunctionBase displacementImposed;
314 
315  nonZero.setFunction (bcNonZero);
316  displacementImposed.setFunction( analyticalDisplacement );
317 
318  //! =================================================================================
319  //! BC for StructuredCube4_test_structuralsolver.mesh
320  //! =================================================================================
321  BCh->addBC ("EdgesIn", 20, Natural, Component, nonZero, compy);
322  BCh->addBC ("EdgesIn", 40, Essential, Component, zero, compy);
323 
324  //! Symmetry BC
325  BCh->addBC ("EdgesIn", 50, EssentialVertices, Component, zero, compxy);
326  BCh->addBC ("EdgesIn", 30, EssentialVertices, Component, zero, compyz);
327  BCh->addBC ("EdgesIn", 80, EssentialVertices, Component, zero, compxz);
328  BCh->addBC ("EdgesIn", 100, EssentialVertices, Full, zero, 3);
329 
330  BCh->addBC ("EdgesIn", 7, Essential, Component, zero, compx);
331  BCh->addBC ("EdgesIn", 3, Essential, Component, zero, compz);
332 
333  //! =================================================================================
334 
335  //! 1. Constructor of the structuralSolver
336  StructuralOperator< RegionMesh<LinearTetra> > solid;
337 
338  //! 2. Setup of the structuralSolver
339  solid.setup (dataStructure,
340  dFESpace,
341  dETFESpace,
342  BCh,
343  parameters->comm);
344 
345  //! 3. Setting data from getPot
346  solid.setDataFromGetPot (dataFile);
347 
348  // Comment this part in order not to have body force in the RHS of the equation
349  bool bodyForce = dataFile ( "solid/physics/bodyForce" , false );
350  if( bodyForce )
351  {
352  solid.setHavingSourceTerm( bodyForce );
353  bodyFunctionPtr_Type bodyTerm( new bodyFunction_Type( f ) );
354 
355  solid.setSourceTerm( bodyTerm );
356  }
357 
358  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
359  {
360  // Setting the vector of fibers functions
361  for( UInt k(1); k <= pointerToVectorOfFamilies->size( ); k++ )
362  {
363  // Setting up the name of the function to define the family
364  std::string family="Family";
365  // adding the number of the family
366  std::string familyNumber;
367  std::ostringstream number;
368  number << ( k );
369  familyNumber = number.str();
370 
371  // Name of the function to create
372  std::string creationString = family + familyNumber;
373  (*pointerToVectorOfFamilies)[ k-1 ].reset( new analyticalFunction_Type() );
374  (*pointerToVectorOfFamilies)[ k-1 ] = setOfFiberFunctions.fiberDefinition( creationString );
375 
376  fiberDirections[ k-1 ].reset( new vector_Type(solid.displacement(), Unique) );
377  }
378 
379  //! 3.b Setting the fibers in the abstract class of Anisotropic materials
380  solid.material()->anisotropicLaw()->setupFiberDirections( pointerToVectorOfFamilies );
381  }
382 
383  //! 4. Building system using TimeAdvance class
384  double timeAdvanceCoefficient = timeAdvance->coefficientSecondDerivative ( 0 ) / (dataStructure->dataTime()->timeStep() * dataStructure->dataTime()->timeStep() );
385  solid.buildSystem (timeAdvanceCoefficient);
386 
387  vectorPtr_Type analyticDispl (new vector_Type ( dFESpace->map() ) );
388  // Interpolating the solution on the mesh
389  dFESpace->interpolate ( static_cast<solidFESpace_Type::function_Type> ( analyticalDisplacement ), *analyticDispl, 0.0 );
390 
391 
392  //! =================================================================================
393  //! Temporal data and initial conditions
394  //! =================================================================================
395 
396  //! 5. Initial data
397  Real dt = dataStructure->dataTime()->timeStep();
398 
399  vectorPtr_Type rhs (new vector_Type ( dFESpace->map() ) );
400  vectorPtr_Type disp (new vector_Type ( dFESpace->map() ) );
401  vectorPtr_Type vel (new vector_Type ( dFESpace->map() ) );
402  vectorPtr_Type acc (new vector_Type ( dFESpace->map() ) );
403 
404  std::vector<vectorPtr_Type> uv0;
405 
406  vectorPtr_Type initialDisplacement (new vector_Type (solid.displacement(), Unique) );
407 
408  if (timeAdvanceMethod == "BDF")
409  {
410  Real tZero = dataStructure->dataTime()->initialTime();
411 
412  for ( UInt previousPass = 0; previousPass < timeAdvance->size() ; previousPass++)
413  {
414  Real previousTimeStep = tZero - previousPass * dt;
415  std::cout << "BDF " << previousTimeStep << "\n";
416  uv0.push_back (disp);
417  }
418  }
419 
420  timeAdvance->setInitialCondition (uv0);
421 
422  timeAdvance->setTimeStep ( dt );
423 
424  timeAdvance->updateRHSContribution ( dt );
425 
426  solid.initialize ( disp );
427 
428  MPI_Barrier (MPI_COMM_WORLD);
429 
430  std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporter;
431  std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporterCheck;
432 
433  std::string const exporterType = dataFile ( "exporter/type", "ensight");
434  std::string const exporterNameFile = dataFile ( "exporter/nameFile", "structure");
435  std::string const exporterCheckName = dataFile ( "exporter/nameCheck", "verifyVectors");
436 
437 #ifdef HAVE_HDF5
438  if (exporterType.compare ("hdf5") == 0)
439  {
440  exporter.reset ( new ExporterHDF5<RegionMesh<LinearTetra> > ( dataFile, exporterNameFile ) );
441  exporterCheck.reset ( new ExporterHDF5<RegionMesh<LinearTetra> > ( dataFile, exporterCheckName ) );
442  }
443  else
444 #endif
445  {
446  if (exporterType.compare ("none") == 0)
447  {
448  exporter.reset ( new ExporterEmpty<RegionMesh<LinearTetra> > ( dataFile, meshPart.meshPartition(), "structure", parameters->comm->MyPID() ) );
449  }
450 
451  else
452  {
453  exporter.reset ( new ExporterEnsight<RegionMesh<LinearTetra> > ( dataFile, meshPart.meshPartition(), "structure", parameters->comm->MyPID() ) );
454  }
455  }
456 
457  exporter->setPostDir ( "./" );
458  exporter->setMeshProcId ( meshPart.meshPartition(), parameters->comm->MyPID() );
459  exporterCheck->setPostDir ( "./" );
460  exporterCheck->setMeshProcId ( meshPart.meshPartition(), parameters->comm->MyPID() );
461 
462  vectorPtr_Type imposedSolution ( new vector_Type ( dFESpace->map() ) );
463  // vectorPtr_Type bodyForceVector ( new vector_Type ( dFESpace->map() ) );
464 
465  vectorPtr_Type solidDisp ( new vector_Type ( dFESpace->map() ) );
466  vectorPtr_Type solidVel ( new vector_Type ( dFESpace->map() ) );
467  vectorPtr_Type solidAcc ( new vector_Type ( dFESpace->map() ) );
468 
469  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "displacement", dFESpace, solidDisp, UInt (0) );
470  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "velocity", dFESpace, solidVel, UInt (0) );
471  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "acceleration", dFESpace, solidAcc, UInt (0) );
472 
473  //exporterCheck->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "forcingTerm", dFESpace, bodyForceVector, UInt (0) );
474  exporterCheck->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "imposedSolution", dFESpace, imposedSolution, UInt (0) );
475 
476  *imposedSolution = *analyticDispl;
477 
478  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
479  {
480  // Adding the fibers vectors
481  // Setting the vector of fibers functions
482  for( UInt k(1); k <= pointerToVectorOfFamilies->size( ); k++ )
483  {
484  // Setting up the name of the function to define the family
485  std::string family="Family-";
486  // adding the number of the family
487  std::string familyNumber;
488  std::ostringstream number;
489  number << ( k );
490  familyNumber = number.str();
491 
492  // Name of the function to create
493  std::string creationString = family + familyNumber;
494  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, creationString, dFESpace, fiberDirections[ k-1 ], UInt (0) );
495 
496  // Extracting the fibers vectors
497  *(fiberDirections[ k-1 ]) = solid.material()->anisotropicLaw()->ithFiberVector( k );
498  }
499  }
500 
501  exporter->postProcess ( 0 );
502  exporterCheck->postProcess ( 0 );
503  std::cout.precision(16);
504 
505  // Saving the analyitical displacement
506  exporterCheck->postProcess ( 1.0 );
507 
508  Real normVect;
509 
510  // Computing the stiffness vector to set as RHS
511  // solid.material()->computeStiffness( *analyticDispl, 0, 1.0, dataStructure,
512  // solid.mapMarkersVolumes(), solid.mapMarkersIndexes(),
513  // solid.displayerPtr() );
514 
515  // vectorPtr_Type stiffnessVector;
516  // stiffnessVector.reset( new vector_Type( dFESpace->map() ) );
517 
518  // *stiffnessVector = *( solid.material()->stiffVector() );
519  //! =============================================================================
520  //! Temporal loop
521  //! =============================================================================
522  // for (Real time = dt; time <= T; time += dt)
523  for (dataStructure->dataTime()->setTime ( dt ) ; dataStructure->dataTime()->canAdvance( ); dataStructure->dataTime()->updateTime( ) )
524  {
525 
526  //! 6. Updating right-hand side
527  *rhs *= 0;
528  timeAdvance->updateRHSContribution ( dt );
529  *rhs += *solid.massMatrix() * timeAdvance->rhsContributionSecondDerivative() / timeAdvanceCoefficient;
530 
531  // Summing the term coming from the imposed solution
532  //*rhs += *stiffnessVector;
533 
534  if( !solid.havingSourceTerm() )
535  {
536  solid.setRightHandSide ( *rhs );
537  }
538  else
539  {
540  solid.updateRightHandSideWithBodyForce( dataStructure->dataTime()->time(), *rhs );
541  }
542 
543  //debug
544  //*bodyForceVector = solid.bodyForce();
545 
546  //! 7. Iterate --> Calling Newton
547  solid.iterate ( BCh );
548 
549  timeAdvance->shiftRight ( solid.displacement() );
550 
551  *solidDisp = solid.displacement();
552  *solidVel = timeAdvance->firstDerivative();
553  *solidAcc = timeAdvance->secondDerivative();
554 
555  exporter->postProcess ( dataStructure->dataTime()->time() );
556  //!--------------------------------------------------------------------------------------------------
557 
558  MPI_Barrier (MPI_COMM_WORLD);
559  }
560  Real errorAbs, error;
561 
562  // Computing the L2 error of the solution
563  errorAbs = dFESpace->l2Error( static_cast<solidFESpace_Type::function_Type> ( analyticalDisplacement ),
564  solid.displacement(), 0.0, &error);
565 
566  std::cout << "L2 relative error: " << error << std::endl;
567  std::cout << "L2 abs error: " << errorAbs << std::endl;
568 }
569 
570 int
571 main ( int argc, char** argv )
572 {
573 
574 #ifdef HAVE_MPI
575  MPI_Init (&argc, &argv);
576  std::shared_ptr<Epetra_MpiComm> Comm (new Epetra_MpiComm ( MPI_COMM_WORLD ) );
577  if ( Comm->MyPID() == 0 )
578  {
579  std::cout << "% using MPI" << std::endl;
580  }
581 #else
582  std::shared_ptr<Epetra_SerialComm> Comm ( new Epetra_SerialComm() );
583  std::cout << "% using serial Version" << std::endl;
584 #endif
585 
586  Structure structure ( argc, argv, Comm );
587  structure.run();
588 
589 #ifdef HAVE_MPI
590  MPI_Finalize();
591 #endif
592  return returnValue ;
593 }
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
Structure(int argc, char **argv, std::shared_ptr< Epetra_Comm > structComm)
std::function< VectorSmall< 3 > Real const &, Real const &, Real const &, Real const &) > bodyFunction_Type
int main(int argc, char **argv)
GetPot(const int argc_, char **argv_, const char *FieldSeparator=0x0)
Definition: GetPot.hpp:507
std::shared_ptr< hdf5Filter_Type > hdf5FilterPtr_Type
FESpace< RegionMesh< LinearTetra >, MapEpetra > solidFESpace_Type
std::shared_ptr< analyticalFunction_Type > analyticalFunctionPtr_Type
void run3d()
run the 3D driven cylinder simulation
ETFESpace< RegionMesh< LinearTetra >, MapEpetra, 3, 3 > solidETFESpace_Type
void run2d()
run the 2D driven cylinder simulation
MeshData - class for handling spatial discretization.
Definition: MeshData.hpp:72
std::shared_ptr< bodyFunction_Type > bodyFunctionPtr_Type
std::shared_ptr< vector_Type > vectorPtr_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
std::set< UInt > parseList(const std::string &list)
std::shared_ptr< solidETFESpace_Type > solidETFESpacePtr_Type
std::function< Real(Real const &, Real const &, Real const &, Real const &, ID const &) > analyticalFunction_Type
std::shared_ptr< TimeAdvance< vector_Type > > timeAdvance_Type
std::vector< vectorPtr_Type > listOfFiberDirections_Type
bool operator()(const char *VarName, bool Default) const
Definition: GetPot.hpp:2008
StructuralOperator< RegionMesh< LinearTetra > >::vector_Type vector_Type
std::shared_ptr< solidFESpace_Type > solidFESpacePtr_Type
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
std::shared_ptr< vectorFiberFunction_Type > vectorFiberFunctionPtr_Type
std::vector< fiberFunctionPtr_Type > vectorFiberFunction_Type