Go to the documentation of this file.
12#ifndef OPM_ELASTICITY_IMPL_HPP
13#define OPM_ELASTICITY_IMPL_HPP
20 template< class Gr idType>
21 template< int components, int funcdim>
23 const Dune::FieldVector<ctype,dim>& point,
24 const Dune::FieldMatrix<ctype,dim,dim>& Jinv)
28 int funcs = funcdim/dim;
30 Dune::FieldMatrix< ctype,funcdim/dim,dim> dNdX;
31 for ( int i=0;i < funcs; i++)
32 Jinv.mv(basis[i].evaluateGradient(point),dNdX[i]);
33 static const int rows = 6+(dim-2)*12;
35 Dune::FieldMatrix< ctype,rows,funcdim/dim> temp;
39#define INDEX(i,j) i+6*j
40 for ( int i=0; i < funcs; ++i) {
42 temp[ INDEX(0,0)][i] = dNdX[i][0];
43 temp[ INDEX(1,1)][i] = dNdX[i][1];
44 temp[ INDEX(2,2)][i] = dNdX[i][2];
47 temp[ INDEX(3,0)][i] = dNdX[i][1];
48 temp[ INDEX(3,1)][i] = dNdX[i][0];
50 temp[ INDEX(4,0)][i] = dNdX[i][2];
51 temp[ INDEX(4,2)][i] = dNdX[i][0];
52 temp[ INDEX(5,1)][i] = dNdX[i][2];
53 temp[ INDEX(5,2)][i] = dNdX[i][1];
55 } else if (dim == 2) {
57#define INDEX(i,j) i+3*j
58 for ( int i=0; i < funcs; ++i) {
60 temp[ INDEX(0,0)][i] = dNdX[i][0];
61 temp[ INDEX(1,1)][i] = dNdX[i][1];
64 temp[ INDEX(2,0)][i] = dNdX[i][1];
65 temp[ INDEX(2,1)][i] = dNdX[i][0];
69 for ( int j=0; j<funcs*dim;++ j)
70 for ( size_t i=0;i<components;++i,++k)
71 B[i][ j] = temp[k%rows][k/rows];
75 template< class Gr idType>
76 template< int funcdims>
78 const Dune::FieldVector<ctype,dim>& point)
83 Dune::FieldMatrix<ctype,funcdims,dim> N;
84 for ( int i = 0; i < funcdims; i++) {
85 Bvector[i] = basis[i].evaluateFunction(point);
89 template < class Gr idType>
90 template< int comp, int funcdim>
92 Dune::FieldMatrix<ctype,funcdim,funcdim>& A,
93 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
94 const Dune::FieldMatrix<ctype,comp,comp>& C,
97 Dune::FieldMatrix<ctype,funcdim,comp> B2;
98 for ( int i=0;i<comp;++i)
99 for ( int j=0; j<funcdim;++ j)
101 A = B.leftmultiplyany(C).leftmultiplyany(B2);
105 template< class Gr idType>
106 template< int comp, int funcdim>
108 const Dune::FieldVector<ctype,funcdim>& v,
109 const Dune::FieldVector<ctype,comp>& eps0,
110 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
111 const Dune::FieldMatrix<ctype,comp,comp>& C)
113 sigma = Dune::FMatrixHelp::mult(C,Dune::FMatrixHelp::mult(B,v)+ eps0);
116Dune::FieldVector<double,3> waveSpeeds( const Dune::FieldMatrix<double,6,6>& C, double phi, double theta, double density)
119 Dune::FieldVector<double, 3> x;
120 x[0] = r*cos(theta)*cos(phi);
121 x[1] = r*sin(theta)*cos(phi);
124 Dune::FieldMatrix<double, 3, 6> D;
135 Dune::FieldMatrix<double, 6, 3> Dt;
136 for ( size_t i=0;i<6;++i)
137 for ( size_t j=0; j<3;++ j)
140 Dune::FieldMatrix<double, 3, 6> T;
141 Dune::FieldMatrix<double, 3, 3> E;
142 Dune::FMatrixHelp::multMatrix(D, C, T);
143 Dune::FMatrixHelp::multMatrix(T, Dt, E);
144 Dune::FieldVector<double, 3> eigenvalues;
145 Dune::FMatrixHelp::eigenValues(E, eigenvalues);
146 std::sort(eigenvalues.begin(), eigenvalues.end(), std::greater<double>());
147 for ( size_t i=0;i<3;++i)
148 eigenvalues[i] = sqrt(eigenvalues[i]/density);
Elasticity helper class. Definition: elasticity.hpp:22
GridType::LeafGridView::ctype ctype A basic number. Definition: elasticity.hpp:28
Singleton handler for the set of LinearShapeFunction. Definition: shapefunctions.hpp:181
static const P1ShapeFunctionSet & instance() Get the only instance of this class. Definition: shapefunctions.hpp:193
int j Definition: elasticity_upscale_impl.hpp:658
Dune::FieldVector< double, 3 > waveSpeeds(const Dune::FieldMatrix< double, 6, 6 > &C, double phi, double theta, double density) Compute the elastic wave velocities. Definition: elasticity_impl.hpp:116
Dune::FieldVector< ctype, comp > eps0 Definition: elasticity_upscale_impl.hpp:446
Definition: ImplicitAssembly.hpp:43
Classes for shape functions. Loosely based on code in dune-grid-howto.
|