LifeV
structure/fem/AssemblyElementalStructure.hpp
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 /*!
28  @file
29  @brief File containing the procedures for the local assembly of the differential operators
30 
31  @author Gianmarco Mengaldo <gianmarco.mengaldo@gmail.com>
32  @author Paolo Tricerri <gianmarco.mengaldo@gmail.com>
33  @mantainer Paolo Tricerri <paolo.tricerri@epfl.ch>
34 
35  All the methods are described in the report StructuralSolver framework in LifeV: Description and Usage
36  */
37 
38 
39 #ifndef _ELEMOPERSTRUCTURE_H_INCLUDED
40 #define _ELEMOPERSTRUCTURE_H_INCLUDED
41 
42 #include <string>
43 #include <iostream>
44 
45 //Trilinos includ
46 #include <Epetra_LAPACK.h>
47 #include <Epetra_BLAS.h>
48 #include <Epetra_SerialDenseMatrix.h>
49 #include <Epetra_SerialDenseVector.h>
50 
51 #include <lifev/core/array/MatrixElemental.hpp>
52 #include <lifev/core/array/VectorElemental.hpp>
53 
54 #include <lifev/core/LifeV.hpp>
55 
56 #include <lifev/core/fem/FESpace.hpp>
57 #include <lifev/core/fem/CurrentFEManifold.hpp>
58 #include <lifev/core/fem/CurrentFE.hpp>
59 #include <lifev/core/fem/DOF.hpp>
60 
61 #include <lifev/core/array/VectorEpetra.hpp>
62 
63 namespace LifeV
64 {
65 //! @name Public typedefs
66 //@{
67 typedef boost::numeric::ublas::matrix<Real> Matrix;
68 typedef boost::numeric::ublas::vector<Real> Vector;
69 typedef boost::numeric::ublas::zero_matrix<Real> ZeroMatrix;
70 //@}
71 
72 /*! /namespace AssemblyElementalStructure
73 
74  This namespace is specially designed to contain the elementary
75  operations (corresponding to differential operators) that build
76  the local contributions to be used in the assembly procedures.
77 
78 */
79 namespace AssemblyElementalStructure
80 {
81 
82 //! Save displacement according to a functor
83 /*!
84  @param FE the finite element space
85  @param fct the functor for the saving
86  @param originVector the vector from which the values are taken
87  @param saveVector the vector from which the values are saved
88 */
89 template<typename FunctorType, typename MeshType, typename MapType>
90 void saveVectorAccordingToFunctor ( const std::shared_ptr<FESpace<MeshType, MapType> > dispFESpace,
91  const FunctorType functor,
92  const VectorEpetra& originVector,
93  const std::shared_ptr<VectorEpetra> statusVector,
94  const std::shared_ptr<VectorEpetra> saveVector,
95  const UInt offset)
96 {
97  UInt dim = dispFESpace->dim();
98  UInt nTotDof = dispFESpace->dof().numTotalDof();
99 
100  // We loop over the local ID on the processors of the originVector
101  for( UInt i(0); i < nTotDof; ++i )
102  {
103  if( originVector.blockMap().LID( static_cast<EpetraInt_Type>(i) ) != -1 ) // The i-th dof is on the processor
104  {
105  // We look at the functor to see if the condition is satified
106  bool variable = functor( i, nTotDof, offset );
107 
108  if ( variable )
109  {
110  // At the first time we insert the value and then we "close" the cell
111  for ( UInt iComp = 0; iComp < dispFESpace->fieldDim(); ++iComp )
112  {
113  Int LIDid = originVector.blockMap().LID (static_cast<EpetraInt_Type>(i + iComp * dim + offset));
114  Int GIDid = originVector.blockMap().GID (LIDid);
115 
116  if( (*statusVector)( GIDid ) != 1.0 )
117  {
118  // Saving the value
119  (*saveVector)( GIDid ) = originVector( GIDid );
120 
121  // Saying it has been saved at that point
122  (*statusVector) ( GIDid ) = 1.0;
123  }
124  }
125  }
126  }
127  }
128 }
129 
130 template<typename FunctorType, typename MeshType, typename MapType>
131 void saveBooleanVectorAccordingToFunctor ( const std::shared_ptr<FESpace<MeshType, MapType> > dispFESpace,
132  const FunctorType& functor,
133  const std::shared_ptr<VectorEpetra> originVector,
134  std::shared_ptr<VectorEpetra> saveVector,
135  const UInt offset)
136 {
137  UInt dim = dispFESpace->dim();
138  UInt nTotDof = dispFESpace->dof().numTotalDof();
139 
140  // We loop over the local ID on the processors of the originVector
141  for( UInt i(0); i < nTotDof; ++i )
142  {
143  if( originVector->blockMap().LID( static_cast<EpetraInt_Type>(i) ) != -1 ) // The i-th dof is on the processor
144  {
145  // We look at the functor to see if the condition is satified
146  bool variable = functor( i, nTotDof, offset );
147 
148  if ( variable )
149  {
150  // At the first time we insert the value and then we "close" the cell
151  for ( UInt iComp = 0; iComp < dispFESpace->fieldDim(); ++iComp )
152  {
153  Int LIDid = originVector->blockMap().LID (static_cast<EpetraInt_Type>(i + iComp * dim + offset));
154  Int GIDid = originVector->blockMap().GID (LIDid);
155 
156  (*saveVector)( GIDid ) = 1.0;
157  }
158  }
159  }
160  }
161 }
162 
163 //! Gradient of the displacement on the local element
164 /*!
165  This function assembles the local tensor of the gradient of the displacement field
166 
167  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
168  @param fe The current finite element
169 */
170 void computeGradientLocalDisplacement (boost::multi_array<Real, 3>& gradientLocalDisplacement, const VectorElemental& uk_loc, const CurrentFE& fe );
171 
172 //! Deformation Gradient on the local element
173 /*!
174  This function assembles the local deformation gradient
175 
176  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
177  @param fe The current finite element
178 */
179 void computeLocalDeformationGradient (const VectorElemental& uk_loc, std::vector<Epetra_SerialDenseMatrix>& tensorF, const CurrentFE& fe );
180 
181 //! Gradient on the local element
182 /*!
183  This function assembles the local gradient
184 
185  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
186  @param fe The current finite element
187 */
188 void computeLocalDeformationGradientWithoutIdentity (const VectorElemental& uk_loc, std::vector<Epetra_SerialDenseMatrix>& tensorF, const CurrentFE& fe );
189 
190 
191 //! METHODS SHARED BETWEEN LINEAR ELASTIC MODEL AND ST.VENANT-KIRCHHOFF MODEL
192 //! These two methods are implemented in AssemblyElemental.cpp.
193 //void stiff_strain( Real coef, MatrixElemental& elmat, const CurrentFE& fe );
194 //void stiff_div( Real coef, MatrixElemental& elmat, const CurrentFE& fe );
195 
196 
197 //! METHODS FOR ST.VENANT-KIRCHHOFF MODEL
198 //! Methods for the Stiffness matrix ( evaluate the RHS or the residuals )
199 
200 //! Elementary first term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
201 /*!
202  This function assembles the local first term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
203 
204  @param coef The constant coefficient of the matrix
205  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
206  @param elmat The elementary matrix of the current volume
207  @param fe The current finite element
208 */
209 void stiff_derdiv ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
210 
211 
212 //! Elementary second term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
213 /*!
214  This function assembles the local second term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
215 
216  @param coef The constant coefficient of the matrix
217  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
218  @param elmat The elementary matrix of the current volume
219  @param fe The current finite element
220 */
221 void stiff_dergradbis ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
222 
223 
224 //! Elementary third term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
225 /*!
226  This function assembles the local third term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
227 
228  @param coef The constant coefficient of the matrix
229  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
230  @param elmat The elementary matrix of the current volume
231  @param fe The current finite element
232 */
233 void stiff_divgrad ( Real coef, const VectorElemental& uk_loc, MatrixElemental& elmat, const CurrentFE& fe );
234 
235 
236 //! Elementary fourth term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
237 /*!
238  This function assembles the local fourth term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
239 
240  @param coef The constant coefficient of the fourth term of the nonlinear stiffness matrix
241  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
242  @param elmat The elementary stiffness matrix of the current volume
243  @param fe The current finite element
244 */
245 void stiff_gradgrad ( Real coef, const VectorElemental& uk_loc, MatrixElemental& elmat, const CurrentFE& fe );
246 
247 
248 //! Elementary fifth term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
249 /*!
250  This function assembles the local fifth term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
251 
252  @param coef The constant coefficient of the matrix
253  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
254  @param elmat The elementary matrix of the current volume
255  @param fe The current finite element
256 */
257 void stiff_dergrad_gradbis ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
258 
259 
260 //! Elementary fifth-2 term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
261 /*!
262  This function assembles the local fifth-2 term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
263 
264  @param coef The constant coefficient of the matrix
265  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
266  @param elmat The elementary matrix of the current volume
267  @param fe The current finite element
268 */
269 void stiff_dergrad_gradbis_Tr ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
270 
271 
272 //! Elementary sixth term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model (see the reference)
273 /*!
274  This function assembles the local sixth term of the nonlinear stiffness matrix for St.Venant-Kirchhoff model.
275 
276  @param coef The constant coefficient of the matrix
277  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
278  @param elmat The elementary matrix of the current volume
279  @param fe The current finite element
280 */
281 void stiff_gradgradTr_gradbis ( Real coef, const VectorElemental& uk_loc, MatrixElemental& elmat, const CurrentFE& fe );
282 
283 
284 
285 
286 //! METHODS FOR THE JACOBIAN MATRIX
287 
288 //! Elementary first term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
289 /*!
290  This function assembles the local first term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
291 
292  @param coef The constant coefficient of the matrix
293  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
294  @param elmat The elementary matrix of the current volume
295  @param fe The current finite element
296 */
297 void stiff_dergrad ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
298 
299 
300 //! Elementary second term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
301 /*!
302  This function assembles the local second term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
303 
304  @param coef The constant coefficient of the matrix
305  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
306  @param elmat The elementary matrix of the current volume
307  @param fe The current finite element
308 */
309 void stiff_divgrad_2 ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
310 
311 
312 //! Elementary third term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
313 /*!
314  This function assembles the local third term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
315 
316  @param coef The constant coefficient of the matrix
317  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
318  @param elmat The elementary matrix of the current volume
319  @param fe The current finite element
320 */
321 void stiff_gradgrad_2 ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
322 
323 
324 //! Elementary fourth term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
325 /*!
326  This function assembles the local fourth term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
327 
328  @param coef The constant coefficient of the matrix
329  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
330  @param elmat The elementary matrix of the current volume
331  @param fe The current finite element
332 */
333 void stiff_dergrad_gradbis_2 ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
334 
335 
336 //! Elementary fifth term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
337 /*!
338  This function assembles the local fifth term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
339 
340  @param coef The constant coefficient of the matrix
341  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
342  @param elmat The elementary matrix of the current volume
343  @param fe The current finite element
344 */
345 void stiff_dergrad_gradbis_Tr_2 ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
346 
347 
348 //! Elementary sixth term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
349 /*!
350  This function assembles the local sixth term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
351 
352  @param coef The constant coefficient of the matrix
353  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
354  @param elmat The elementary matrix of the current volume
355  @param fe The current finite element
356 */
357 void stiff_gradgradTr_gradbis_2 ( Real coef, const boost::multi_array<Real, 3>& gradientLocalDisplacement, MatrixElemental& elmat, const CurrentFE& fe );
358 
359 
360 //! Elementary seventh term of the Jacobian matrix for the nonlinear stiffness matrix of the St.Venant-Kirchhoff model (see the reference)
361 /*!
362  This function assembles the local seventh term of the Jacobian matrix of the nonlinear stiffness matrix of the St.Venant-Kirchhoff model.
363 
364  @param coef The constant coefficient of the matrix
365  @param uk_loc The local displacement (remark: the nonlinear matrix depends on current displacement)
366  @param elmat The elementary matrix of the current volume
367  @param fe The current finite element
368 */
369 void stiff_gradgradTr_gradbis_3 ( Real coef, const VectorElemental& uk_loc, MatrixElemental& elmat, const CurrentFE& fe );
370 
371 
372 
373 //! METHODS SHARED BETWEEN NEO-HOOKEAN AND EXPONENTIAL MODELS
374 //! The volumetric part for Neo-Hookean and Exponential models is the same
375 //! Methods for the volumetric part of the stiffness vector
376 
377 //! Elementary volumetric term of the nonlinear stiffness vector of the Neo-Hookean and Exponential models (see the reference)
378 /*!
379  This function assembles the volumetric term of the nonlinear stiffness vector of the Neo-Hookean and Exponential models.
380 
381  @param coef The constant coefficient of the volumetric term of the nonlinear stiffness vector
382  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
383  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
384  @param elvec The elementary vector of the current volume
385  @param fe The current finite element
386 */
387 void source_Pvol (Real coef, const boost::multi_array<Real, 3>& CofFk, const std::vector<Real>& Jk, VectorElemental& elvec, const CurrentFE& fe);
388 
389 //! Methods for the volumetric part of the Jacobian matrix
390 
391 //! Elementary first volumetric term of the nonlinear Jacobian matrix of the Neo-Hookean and Exponential models (see the reference)
392 /*!
393  This function assembles the local first volumetric term of the Jacobian matrix of the nonlinear volumetric stiffness vector of the Neo-Hookean and Exponential models.
394 
395  @param coef The constant coefficient of the matrix
396  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
397  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
398  @param elvec The elementary matrix of the current volume
399  @param fe The current finite element
400 */
401 void stiff_Jac_Pvol_1term ( Real coef, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
402 
403 
404 //! Elementary second volumetric term of the nonlinear Jacobian matrix of the Neo-Hookean and Exponential models (see the reference)
405 /*!
406  This function assembles the local second volumetric term of the Jacobian matrix of the nonlinear volumetric stiffness vector of the Neo-Hookean and Exponential models.
407 
408  @param coef The constant coefficient of the matrix
409  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
410  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
411  @param elmat The elementary matrix of the current volume
412  @param fe The current finite element
413 */
414 void stiff_Jac_Pvol_2term ( Real coef, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
415 
416 //! METHODS FOR NEO-HOOKEAN MODEL
417 //! Methods for the isochoric part of the stiffness vector
418 
419 //! Elementary nonlinear isochoric stiffness vector for Neo-Hookean model (see the reference)
420 /*!
421  This function assembles the nonlinear isochoric part of the stiffness vector for Neo-Hookean model.
422 
423  @param coef The coefficient of the vector
424  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
425  @param Fk The deformation gradient that depends on the local displacement uk_loc
426  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
427  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
428  @param elvec The elementary vector of the current volume
429  @param fe The current finite element
430 */
431 void source_P1iso_NH (Real coef, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, VectorElemental& elvec, const CurrentFE& fe);
432 
433 //! METHODS FOR TENSORIAL CALCULUS
434 //! In this part of the namespace, the methods to perform basics operations on tensors are defined.
435 //! The operations are: tensorial products, computations of invariants
436 /*!
437  This function computes the invariants of the right Cauchy Green tensor and the cofactor of F
438 
439  @param invariants vector of invariants of C
440  @param tensorF deformation gradient tensor
441 */
442 void computeInvariantsRightCauchyGreenTensor (std::vector<LifeV::Real>& invariants,
443  const Epetra_SerialDenseMatrix& tensorF,
444  Epetra_SerialDenseMatrix& cofactorF);
445 
446 void computeInvariantsRightCauchyGreenTensor (std::vector<LifeV::Real>& invariants,
447  const Epetra_SerialDenseMatrix& tensorF);
448 
449 
450 /*!
451  This function computes the Cauchy stress tensor given the detF, first Piola Kirchhoff and tensorF
452 
453  @param cauchy Cauchy stress tensor
454  @param firstPiola first Piola-Kirchhoff tensor
455  @param invariants vector of invariants of C
456  @param tensorF deformation gradient tensor
457 */
458 void computeCauchyStressTensor (Epetra_SerialDenseMatrix& cauchy,
459  Epetra_SerialDenseMatrix& firstPiola,
460  LifeV::Real det,
461  Epetra_SerialDenseMatrix& tensorF);
462 
463 /*!
464  This function computes the eigenvalues of \sigma
465  @param cauchy Cauchy stress tensor
466  @param eigenvalues vector of principal tensions
467 */
468 void computeEigenvalues (const Epetra_SerialDenseMatrix& cauchy,
469  std::vector<LifeV::Real>& eigenvaluesR,
470  std::vector<LifeV::Real>& eigenvaluesI);
471 
472 //! Methods for the isochoric part of the Jacobian matrix
473 
474 //! Elementary first nonlinear isochoric Jacobian matrix for Neo-Hookean model (see the reference)
475 /*!
476  This function assembles the local first nonlinear isochooric Jacobian matrix for Neo-Hookean model.
477 
478  @param coef The coefficient of the matrix
479  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
480  @param Fk The deformation gradient that depends on the local displacement uk_loc
481  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
482  @param elmat The elementary matrix of the current volume
483  @param fe The current finite element
484 */
485 void stiff_Jac_P1iso_NH_1term ( Real coef, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
486 
487 //! Elementary second nonlinear isochoric Jacobian matrix for Neo-Hookean model (see the reference)
488 /*!
489  This function assembles the local second nonlinear isochooric Jacobian matrix for Neo-Hookean model.
490 
491  @param coef The coefficient of the matrix
492  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
493  @param Fk The deformation gradient that depends on the local displacement uk_loc
494  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
495  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
496  @param elmat The elementary matrix of the current volume
497  @param fe The current finite element
498 */
499 void stiff_Jac_P1iso_NH_2term ( Real coef, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
500 
501 
502 //! Elementary third nonlinear isochoric Jacobian matrix for Neo-Hookean model (see the reference)
503 /*!
504  This function assembles the local third nonlinear isochooric Jacobian matrix for Neo-Hookean model.
505 
506  @param coef The coefficient of the matrix
507  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
508  @param elmat The elementary matrix of the current volume
509  @param fe The current finite element
510 */
511 void stiff_Jac_P1iso_NH_3term ( Real coef, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
512 
513 
514 //! Elementary fourth nonlinear isochoric Jacobian matrix for Neo-Hookean model (see the reference)
515 /*!
516  This function assembles the local fourth nonlinear isochooric Jacobian matrix for Neo-Hookean model.
517 
518  @param coef The coefficient of the matrix
519  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
520  @param Fk The deformation gradient that depends on the local displacement uk_loc
521  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
522  @param elmat The elementary matrix of the current volume
523  @param fe The current finite element
524 */
525 void stiff_Jac_P1iso_NH_4term ( Real coef, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
526 
527 
528 //! Elementary fifth nonlinear isochoric Jacobian matrix for Neo-Hookean model (see the reference)
529 /*!
530  This function assembles the local fifth nonlinear isochooric Jacobian matrix for Neo-Hookean model.
531 
532  @param coef The coefficient of the matrix
533  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
534  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
535  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
536  @param elmat The elementary matrix of the current volume
537  @param fe The current finite element
538 */
539 void stiff_Jac_P1iso_NH_5term ( Real coef, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
540 
541 //! METHODS FOR EXPONENTIAL MODEL
542 //! Methods for the isochoric part of the stiffness vector
543 
544 //! Elementary nonlinear isochoric stiffness vector for Exponential model (see the reference)
545 /*!
546  This function assembles the local nonlinear isochoric stiffness vector for Exponential model.
547 
548  @param coef The pre-exponential coefficient
549  @param coefExp The expoenential coefficient
550  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
551  @param Fk The deformation gradient that depends on the local displacement uk_loc
552  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
553  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
554  @param elvec The elementary vector of the current volume
555  @param fe The current finite element
556 */
557 void source_P1iso_Exp ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, VectorElemental& elvec, const CurrentFE& fe );
558 
559 //! Methods for the isochoric part of the Jacobian matrix
560 
561 //! Elementary first nonlinear isochoric Jacobian matrix for Exponential model (see the reference)
562 /*!
563  This function assembles the local first nonlinear isochoric Jacobian matrix for Exponential model.
564 
565  @param coef The pre-exponential coefficient
566  @param coefExp The expoenential coefficient
567  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
568  @param Fk The deformation gradient that depends on the local displacement uk_loc
569  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
570  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
571  @param elmat The elementary matrix of the current volume
572  @param fe The current finite element
573 */
574 void stiff_Jac_P1iso_Exp_1term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
575 
576 
577 //! Elementary second nonlinear isochoric Jacobian matrix for Exponential model (see the reference)
578 /*!
579  This function assembles the local second nonlinear isochoric Jacobian matrix for Exponential model.
580 
581  @param coef The pre-exponential coefficient
582  @param coefExp The expoenential coefficient
583  @param Fk The deformation gradient that depends on the local displacement uk_loc
584  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
585  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
586  @param elmat The elementary matrix of the current volume
587  @param fe The current finite element
588 */
589 void stiff_Jac_P1iso_Exp_2term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
590 
591 
592 //! Elementary third nonlinear isochoric Jacobian matrix for Exponential model (see the reference)
593 /*!
594  This function assembles the local third nonlinear isochoric Jacobian matrix for Exponential model.
595 
596  @param coef The pre-exponential coefficient
597  @param coefExp The expoenential coefficient
598  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
599  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
600  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
601  @param elmat The elementary matrix of the current volume
602  @param fe The current finite element
603 */
604 void stiff_Jac_P1iso_Exp_3term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
605 
606 
607 //! Elementary fourth nonlinear isochoric Jacobian matrix for Exponential model (see the reference)
608 /*!
609  This function assembles the local fourth nonlinear isochoric Jacobian matrix for Exponential model.
610 
611  @param coef The pre-exponential coefficient
612  @param coefExp The expoenential coefficient
613  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
614  @param Fk The deformation gradient that depends on the local displacement uk_loc
615  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
616  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
617  @param elmat The elementary matrix of the current volume
618  @param fe The current finite element
619 */
620 void stiff_Jac_P1iso_Exp_4term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
621 
622 //! Elementary fifth nonlinear isochoric Jacobian matrix for Exponential model (see the reference)
623 /*!
624  This function assembles the local fifth nonlinear isochoric Jacobian matrix for Exponential model.
625 
626  @param coef The pre-exponential coefficient
627  @param coefExp The expoenential coefficient
628  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
629  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
630  @param elmat The elementary matrix of the current volume
631  @param fe The current finite element
632 */
633 void stiff_Jac_P1iso_Exp_5term ( Real coef, Real coefExp, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
634 
635 //! Elementary sixth nonlinear isochoric Jacobian matrix for Exponential model (see the reference)
636 /*!
637  This function assembles the local sixth nonlinear isochoric Jacobian matrix for Exponential model.
638 
639  @param coef The pre-exponential coefficient
640  @param coefExp The expoenential coefficient
641  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
642  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
643  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
644  @param elmat The elementary matrix of the current volume
645  @param fe The current finite element
646 */
647 void stiff_Jac_P1iso_Exp_6term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
648 
649 //! METHOD FOR THE ST. VENANT-KIRCHHOFF LAW IN THE ISOCHORIC AND VOLUMETRIC DECOUPLED VERSION
650 //! This law uses, for the volumetric part, the one used by the exponential and neohookean model.
651 //! The isochoric part is, of course, completely new and, consequently its jacobian.
652 
653 //! Elementary nonlinear isochoric stiffness vector for St. Venant-Kirchhoff Penalized model
654 /*!
655  This function assembles the local nonlinear isochoric stiffness vector for St. Venant-Kirchhoff Penalized model
656 
657  @param lambda first Lame Coefficient
658  @param mu shear modulus
659  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
660  @param Fk The deformation gradient that depends on the local displacement uk_loc
661  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
662  @param elvec The elementary vector of the current volume
663  @param fe The current finite element
664 */
665 void source_P1iso_VKPenalized ( Real lambda, Real mu, const boost::multi_array<Real, 3 >& FkMinusTransposed, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Ic_isok, const std::vector<Real>& Ic_k, const std::vector<Real>& Jack_k, VectorElemental& elvec, const CurrentFE& fe );
666 
667 //! Elementary nonlinear isochoric stiffness vector for St. Venant-Kirchhoff Penalized model
668 /*!
669  This function assembles the local nonlinear isochoric stiffness vector for St. Venant-Kirchhoff Penalized model
670 
671  @param mu shear modulus
672  @param Fk^-T
673  @param FkCk The product FC
674  @param Jk The determinant of F
675  @param Ic_Squared The trace of C^2
676  @param elvec The elementary vector of the current volume
677  @param fe The current finite element
678 */
679 void source_P2iso_VKPenalized ( Real mu, const boost::multi_array<Real, 3 >& FkMinusTransposed, const boost::multi_array<Real, 3 >& FkCk, const std::vector<Real>& Ic_Squared, const std::vector<Real>& Jk, VectorElemental& elvec, const CurrentFE& fe );
680 
681 
682 //! Methdos for the Jacobian of the St. Venant-Kirchhoff Penalized law.
683 //! Elementary first nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
684 /*!
685  This function assembles the local first nonlinear isochoric Jacobian matrix for VK-Penalized model.
686 
687  @param lambda the first Lame coefficient
688  @param FkMinusTransposed The matrix F^-T
689  @param Fk The deformation gradient that depends on the local displacement uk_loc
690  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
691  @param Ic_k The trace of C
692  @param elmat The elementary matrix of the current volume
693  @param fe The current finite element
694 */
695 void stiff_Jac_P1iso_VKPenalized_0term ( Real lambda, Real mu, const boost::multi_array<Real, 3 >& FkMinusTransposed, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_k, const std::vector<Real>& IcIso_k, MatrixElemental& elmat, const CurrentFE& fe );
696 
697 //! Elementary first nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
698 /*!
699  This function assembles the local first nonlinear isochoric Jacobian matrix for VK-Penalized model.
700 
701  @param lambda the first Lame coefficient
702  @param FkMinusTransposed The matrix F^-T
703  @param Fk The deformation gradient that depends on the local displacement uk_loc
704  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
705  @param Ic_k The trace of C
706  @param elmat The elementary matrix of the current volume
707  @param fe The current finite element
708 */
709 void stiff_Jac_P1iso_VKPenalized_1term ( Real coeff, const boost::multi_array<Real, 3 >& FkMinusTransposed, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_k, MatrixElemental& elmat, const CurrentFE& fe );
710 
711 //! Elementary third nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
712 /*!
713  This function assembles the local third nonlinear isochoric Jacobian matrix for VK-Penalized model.
714 
715  @param coef lambda / 2 where lambda is the first Lame constant
716  @param FKMinusTransposed The matrix F^-T
717  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
718  @param Ic_k The first invariant of the right Cauchy-Green tensor C
719  @param elmat The elementary matrix of the current volume
720  @param fe The current finite element
721 */
722 void stiff_Jac_P1iso_VKPenalized_2term ( Real coef, const boost::multi_array<Real, 3 >& FkMinusTransposed, const std::vector<Real>& Jk, const std::vector<Real>& Ic_k, MatrixElemental& elmat, const CurrentFE& fe );
723 
724 
725 //! Elementary second nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
726 /*!
727  This function assembles the local second nonlinear isochoric Jacobian matrix for VK-Penalized model.
728 
729  @param coef (lambda/2) where lambda is the first Lame coefficient
730  @param Fk The deformation gradient that depends on the local displacement uk_loc
731  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
732  @param elmat The elementary matrix of the current volume
733  @param fe The current finite element
734 */
735 void stiff_Jac_P1iso_VKPenalized_3term ( Real coef, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
736 
737 //! Elementary fourth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
738 /*!
739  This function assembles the local fourth nonlinear isochoric Jacobian matrix for VK-Penalized model.
740 
741  @param coef (lambda/2) where lamba is the first Lame coefficient
742  @param FkMinusTransposed The tensor F^-T
743  @param Fk The deformation gradient that depends on the local displacement uk_loc
744  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
745  @param Ic_k The first invariant of the right Cauchy-Green tensor C
746  @param elmat The elementary matrix of the current volume
747  @param fe The current finite element
748 */
749 void stiff_Jac_P1iso_VKPenalized_4term ( Real coef, const boost::multi_array<Real, 3 >& FkMinusTransposed, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_k, MatrixElemental& elmat, const CurrentFE& fe );
750 
751 //! Elementary fifth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
752 /*!
753  This function assembles the local fifth nonlinear isochoric Jacobian matrix for VK-Penalized model.
754 
755  @param coef lambda first lame coefficient
756  @param coef2 mu shear modulus
757  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
758  @param elmat The elementary matrix of the current volume
759  @param fe The current finite element
760 */
761 void stiff_Jac_P1iso_VKPenalized_5term ( Real coef, Real secondCoef, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
762 
763 //! Elementary sixth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
764 /*!
765  This function assembles the local sixth nonlinear isochoric Jacobian matrix for VK-Penalized model.
766 
767  @param coef lambda first Lame coefficient
768  @param coef2 mu sheat modulus
769  @param Ic_k The first invariant of the right Cauchy-Green tensor C
770  @param Fk The deformation gradient that depends on the local displacement uk_loc
771  @param elmat The elementary matrix of the current volume
772  @param fe The current finite element
773 */
774 void stiff_Jac_P1iso_VKPenalized_6term ( Real coef, Real secondCoef, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, const boost::multi_array<Real, 3 >& Fk, const boost::multi_array<Real, 3 >& FkMinusTransposed, MatrixElemental& elmat, const CurrentFE& fe );
775 
776 
777 //! Elementary seventh nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
778 /*!
779  This function assembles the local seventh nonlinear isochoric Jacobian matrix for VK-Penalized model.
780 
781  @param coef lambda the first Lame coefficient
782  @param coef2 mu the shear modulus
783  @param FkMinusTransposed The tensor F^-T
784  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
785  @param Ic_k The first invariant of the right Cauchy-Green tensor C
786  @param elmat The elementary matrix of the current volume
787  @param fe The current finite element
788 */
789 void stiff_Jac_P1iso_VKPenalized_7term ( Real coef, Real secondCoef, const boost::multi_array<Real, 3 >& FkMinusTransposed, const std::vector<Real>& Ic_isok, const std::vector<Real>& Ic_k, const std::vector<Real>& Jk, MatrixElemental& elmat, const CurrentFE& fe );
790 
791 //! Elementary seventh nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
792 /*!
793  This function assembles the local seventh nonlinear isochoric Jacobian matrix for VK-Penalized model.
794 
795  @param coef lambda the first Lame coefficient
796  @param coef2 mu the shear modulus
797  @param FkMinusTransposed The tensor F^-T
798  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
799  @param Ic_k The first invariant of the right Cauchy-Green tensor C
800  @param elmat The elementary matrix of the current volume
801  @param fe The current finite element
802 */
803 void stiff_Jac_P1iso_VKPenalized_8term ( Real coef, const std::vector<Real>& Jack_k, const boost::multi_array<Real, 3 >& FkMinusTransposed, const boost::multi_array<Real, 3 >& FkCk, MatrixElemental& elmat, const CurrentFE& fe );
804 
805 
806 //! Elementary seventh nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
807 /*!
808  This function assembles the local seventh nonlinear isochoric Jacobian matrix for VK-Penalized model.
809 
810  @param coef lambda the first Lame coefficient
811  @param coef2 mu the shear modulus
812  @param FkMinusTransposed The tensor F^-T
813  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
814  @param Ic_k The first invariant of the right Cauchy-Green tensor C
815  @param elmat The elementary matrix of the current volume
816  @param fe The current finite element
817 */
818 void stiff_Jac_P1iso_VKPenalized_9term ( Real coef, const std::vector<Real>& Jack_k, const std::vector<Real>& Ic_kSquared, const boost::multi_array<Real, 3 >& FkMinusTransposed, MatrixElemental& elmat, const CurrentFE& fe );
819 
820 
821 //! Elementary eigth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
822 /*!
823  This function assembles the local eigth nonlinear isochoric Jacobian matrix for VK-Penalized model.
824 
825  @param coef mu the shear modulus
826  @param Ck The right Cauchy-Green tensor C
827  @param elmat The elementary matrix of the current volume
828  @param fe The current finite element
829 */
830 void stiff_Jac_P1iso_VKPenalized_10term ( Real coef, const std::vector<Real>& Jack_k, const boost::multi_array<Real, 3 >& Ck, MatrixElemental& elmat, const CurrentFE& fe );
831 
832 //! Elementary sixth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
833 /*!
834  This function assembles the local sixth nonlinear isochoric Jacobian matrix for VK-Penalized model.
835 
836  @param coef mu the shear modulus
837  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
838  @param Fk The deformation tensor
839  @param elmat The elementary matrix of the current volume
840  @param fe The current finite element
841 */
842 void stiff_Jac_P1iso_VKPenalized_11term ( Real coef, const std::vector<Real>& Jk, const boost::multi_array<Real, 3 >& Fk, MatrixElemental& elmat, const CurrentFE& fe );
843 
844 //! Elementary tenth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
845 /*!
846  This function assembles the local tenth nonlinear isochoric Jacobian matrix for VK-Penalized model.
847 
848  @param coef mu the shear modulus
849  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
850  @param Fk The deformation tensor
851  @param elmat The elementary matrix of the current volume
852  @param fe The current finite element
853 */
854 void stiff_Jac_P1iso_VKPenalized_12term ( Real coef, const std::vector<Real>& Jk, const boost::multi_array<Real, 3 >& Fk, MatrixElemental& elmat, const CurrentFE& fe );
855 
856 //! Elementary eleventh nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
857 /*!
858  This function assembles the local eleventh nonlinear isochoric Jacobian matrix for VK-Penalized model.
859 
860  @param coef mu the shear modulus
861  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
862  @param Ic_isok The first invariant of the right Cauchy-Green tensor C^2
863  @param FkMinusTransposed The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
864  @param elmat The elementary matrix of the current volume
865  @param fe The current finite element
866 */
867 void stiff_Jac_P1iso_VKPenalized_13term ( Real coef, const std::vector<Real>& Jk, const std::vector<Real>& Ic_kSquared, const boost::multi_array<Real, 3 >& FkMinusTransposed, MatrixElemental& elmat, const CurrentFE& fe );
868 
869 //! Elementary twelveth nonlinear isochoric Jacobian matrix for VK-Penalized model (see the reference)
870 /*!
871  This function assembles the local twelveth nonlinear isochoric Jacobian matrix for VK-Penalized model.
872 
873  @param coef mu the shear modulus
874  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
875  @param Ck The first invariant of the right Cauchy-Green tensor C^2
876  @param Fk The deformation gradient matrix
877  @param FkMinusTransposed The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
878  @param elmat The elementary matrix of the current volume
879  @param fe The current finite element
880 */
881 void stiff_Jac_P1iso_VKPenalized_14term ( Real coef, const std::vector<Real>& Jk, const boost::multi_array<Real, 3 >& FkCk, const boost::multi_array<Real, 3 >& FkMinusTransposed, MatrixElemental& elmat, const CurrentFE& fe );
882 
883 //! Methods for second order exponential law
884 //! Methods for the first Piola-Kirchhoff tensor
885 void source_P1iso_SecondOrderExponential ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& trCisok, VectorElemental& elvec, const CurrentFE& fe );
886 
887 //! Methods for the Jacobian matrix
888 //! Elementary first nonlinear isochoric Jacobian matrix for Second Order Exponential model (see the reference)
889 /*!
890  This function assembles the local first nonlinear isochoric Jacobian matrix for Exponential model.
891 
892  @param coef The pre-exponential coefficient
893  @param coefExp The expoenential coefficient
894  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
895  @param Fk The deformation gradient that depends on the local displacement uk_loc
896  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
897  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
898  @param elmat The elementary matrix of the current volume
899  @param fe The current finite element
900 */
901 void stiff_Jac_P1iso_SecondOrderExp_1term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
902 
903 
904 //! Elementary second nonlinear isochoric Jacobian matrix for Second Order Exponential model (see the reference)
905 /*!
906  This function assembles the local second nonlinear isochoric Jacobian matrix for Exponential model.
907 
908  @param coef The pre-exponential coefficient
909  @param coefExp The expoenential coefficient
910  @param Fk The deformation gradient that depends on the local displacement uk_loc
911  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
912  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
913  @param elmat The elementary matrix of the current volume
914  @param fe The current finite element
915 */
916 void stiff_Jac_P1iso_SecondOrderExp_2term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
917 
918 //! Elementary third nonlinear isochoric Jacobian matrix for Second Order Exponential model (see the reference)
919 /*!
920  This function assembles the local third nonlinear isochoric Jacobian matrix for Second Order Exponential model.
921 
922  @param coef The pre-exponential coefficient
923  @param coefExp The expoenential coefficient
924  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
925  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
926  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
927  @param elmat The elementary matrix of the current volume
928  @param fe The current finite element
929 */
930 void stiff_Jac_P1iso_SecondOrderExp_3term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
931 
932 //! Elementary fourth nonlinear isochoric Jacobian matrix for Second Order Exponential model (see the reference)
933 /*!
934  This function assembles the local fourth nonlinear isochoric Jacobian matrix for Second Order Exponential model.
935 
936  @param coef The pre-exponential coefficient
937  @param coefExp The expoenential coefficient
938  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
939  @param Fk The deformation gradient that depends on the local displacement uk_loc
940  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
941  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
942  @param Ic_k The first invariant of the isochoric part of the right Cauchy-Green tensor C
943  @param elmat The elementary matrix of the current volume
944  @param fe The current finite element
945 */
946 void stiff_Jac_P1iso_SecondOrderExp_4term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const boost::multi_array<Real, 3 >& Fk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, const std::vector<Real>& Ic_k, MatrixElemental& elmat, const CurrentFE& fe );
947 
948 //! Elementary fifth nonlinear isochoric Jacobian matrix for Second Order Exponential model (see the reference)
949 /*!
950  This function assembles the local fifth nonlinear isochoric Jacobian matrix for Second Order Exponential model.
951 
952  @param coef The pre-exponential coefficient
953  @param coefExp The expoenential coefficient
954  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
955  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
956  @param elmat The elementary matrix of the current volume
957  @param fe The current finite element
958 */
959 void stiff_Jac_P1iso_SecondOrderExp_5term ( Real coef, Real coefExp, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
960 
961 //! Elementary sixth nonlinear isochoric Jacobian matrix for Second Order Exponential model (see the reference)
962 /*!
963  This function assembles the local sixth nonlinear isochoric Jacobian matrix for Second Order Exponential model.
964 
965  @param coef The pre-exponential coefficient
966  @param coefExp The expoenential coefficient
967  @param CofFk The cofactor of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
968  @param Jk The determinant of the deformation gradient F that depends on the local displacement uk_loc (remark: the nonlinear vector depends on current displacement)
969  @param Ic_isok The first invariant of the isochoric part of the right Cauchy-Green tensor C
970  @param elmat The elementary matrix of the current volume
971  @param fe The current finite element
972 */
973 void stiff_Jac_P1iso_SecondOrderExp_6term ( Real coef, Real coefExp, const boost::multi_array<Real, 3 >& CofFk, const std::vector<Real>& Jk, const std::vector<Real>& Ic_isok, MatrixElemental& elmat, const CurrentFE& fe );
974 
975 } //! End namespace AssemblyElementalStructure
976 
977 } //! End namespace LifeV
978 #endif
VectorEpetra - The Epetra Vector format Wrapper.
void updateInverseJacobian(const UInt &iQuadPt)
double Real
Generic real data.
Definition: LifeV.hpp:175
CurrentFE - A primordial class for the assembly of the local matrices/vectors retaining the values on...
Definition: CurrentFE.hpp:243
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191