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