Newer
Older
#ifndef AVERAGE_INTERFACE_HH
#define AVERAGE_INTERFACE_HH
#include <dune/common/fmatrix.hh>
#include <dune/disc/shapefunctions/lagrangeshapefunctions.hh>
#include "../../contact/src/dgindexset.hh"
#include "../../common/crossproduct.hh"
#include "svd.hh"
#include "../linearsolver.hh"

Oliver Sander
committed
// Given a resultant force and torque (from a rod problem), this method computes the corresponding
// Neumann data for a 3d elasticity problem.
template <class GridType>
void computeAveragePressure(const Dune::FieldVector<double,GridType::dimension>& resultantForce,
const Dune::FieldVector<double,GridType::dimension>& resultantTorque,
const BoundaryPatch<GridType>& interface,
const Configuration& crossSection,
Dune::BlockVector<Dune::FieldVector<double, GridType::dimension> >& pressure)
{
const GridType& grid = interface.getGrid();
const int level = interface.level();
const typename GridType::Traits::LevelIndexSet& indexSet = grid.levelIndexSet(level);
const int dim = GridType::dimension;
typedef typename GridType::ctype ctype;
typedef typename GridType::template Codim<dim>::LevelIterator VertexIterator;
// set up output array
DGIndexSet<GridType> dgIndexSet(grid,level);
dgIndexSet.setup(grid,level);
pressure.resize(dgIndexSet.size());
typename GridType::template Codim<0>::LevelIterator eIt = indexSet.template begin<0,Dune::All_Partition>();
typename GridType::template Codim<0>::LevelIterator eEndIt = indexSet.template end<0,Dune::All_Partition>();
for (; eIt!=eEndIt; ++eIt) {
typename GridType::template Codim<0>::Entity::LevelIntersectionIterator nIt = eIt->ilevelbegin();
typename GridType::template Codim<0>::Entity::LevelIntersectionIterator nEndIt = eIt->ilevelend();
for (; nIt!=nEndIt; ++nIt) {
if (!interface.contains(*eIt,nIt))
continue;
const Dune::LagrangeShapeFunctionSet<double, double, dim-1>& baseSet
= Dune::LagrangeShapeFunctions<double, double, dim-1>::general(nIt.intersectionGlobal().type(),1);
if (baseSet.size() != 4)
DUNE_THROW(Dune::NotImplemented, "average interface only for quad faces");
// four rows because a face may have no more than four vertices
Dune::FieldVector<double,4> mu(0);
Dune::FieldVector<double,3> mu_tilde[4][3];
for (int i=0; i<4; i++)
for (int j=0; j<3; j++)
mu_tilde[i][j] = 0;
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
for (int i=0; i<nIt.intersectionGlobal().corners(); i++) {
const Dune::QuadratureRule<double, dim-1>& quad
= Dune::QuadratureRules<double, dim-1>::rule(nIt.intersectionGlobal().type(), dim-1);
for (size_t qp=0; qp<quad.size(); qp++) {
// Local position of the quadrature point
const Dune::FieldVector<double,dim-1>& quadPos = quad[qp].position();
const double integrationElement = nIt.intersectionGlobal().integrationElement(quadPos);
// \mu_i = \int_t \varphi_i \ds
mu[i] += quad[qp].weight() * integrationElement * baseSet[i].evaluateFunction(0,quadPos);
// \tilde{\mu}_i^j = \int_t \varphi_i \times (x - x_0) \ds
Dune::FieldVector<double,dim> worldPos = nIt.intersectionGlobal().global(quadPos);
for (int j=0; j<dim; j++) {
// Vector-valued basis function
Dune::FieldVector<double,dim> phi_i(0);
phi_i[j] = baseSet[i].evaluateFunction(0,quadPos);
mu_tilde[i][j].axpy(quad[qp].weight() * integrationElement,
crossProduct(worldPos-crossSection.r, phi_i));
}
}
}
// std::cout << "tilde{mu}\n" << std::endl;
// for (int i=0; i<4; i++)
// for (int j=0; j<3; j++)
// std::cout << "i: " << i << ", j: " << j << ", " << mu_tilde[i][j] << std::endl;
// Set up matrix
Dune::FieldMatrix<double, 6, 12> matrix(0);
for (int i=0; i<4; i++)
for (int j=0; j<3; j++)
matrix[j][i*3+j] = mu[i];
for (int i=0; i<4; i++)
for (int j=0; j<3; j++)
for (int k=0; k<3; k++)
matrix[3+k][3*i+j] = mu_tilde[i][j][k];
Dune::FieldVector<double,12> u;
Dune::FieldVector<double,6> b;
for (int i=0; i<3; i++) {
b[i] = resultantForce[i];
b[i+3] = resultantTorque[i];
}
// std::cout << b << std::endl;
// std::cout << matrix << std::endl;
//matrix.solve(u,b);
linearSolver(matrix, u, b);
//std::cout << u << std::endl;
for (int i=0; i<3; i++) {
pressure[dgIndexSet(*eIt, nIt.numberInSelf())][i] = u[i];
pressure[dgIndexSet(*eIt, nIt.numberInSelf())+1][i] = u[i+3];
pressure[dgIndexSet(*eIt, nIt.numberInSelf())+2][i] = u[i+6];
pressure[dgIndexSet(*eIt, nIt.numberInSelf())+3][i] = u[i+9];
}

Oliver Sander
committed
// /////////////////////////////////////////////////////////////////////////////////////
// Compute the overall force and torque to see whether the preceding code is correct
// /////////////////////////////////////////////////////////////////////////////////////
Dune::FieldVector<double,3> outputForce(0), outputTorque(0);
eIt = indexSet.template begin<0,Dune::All_Partition>();
eEndIt = indexSet.template end<0,Dune::All_Partition>();

Oliver Sander
committed
for (; eIt!=eEndIt; ++eIt) {
typename GridType::template Codim<0>::Entity::LevelIntersectionIterator nIt = eIt->ilevelbegin();
typename GridType::template Codim<0>::Entity::LevelIntersectionIterator nEndIt = eIt->ilevelend();
for (; nIt!=nEndIt; ++nIt) {
if (!interface.contains(*eIt,nIt))
continue;
const Dune::LagrangeShapeFunctionSet<double, double, dim-1>& baseSet
= Dune::LagrangeShapeFunctions<double, double, dim-1>::general(nIt.intersectionGlobal().type(),1);

Oliver Sander
committed
const Dune::QuadratureRule<double, dim-1>& quad
= Dune::QuadratureRules<double, dim-1>::rule(nIt.intersectionGlobal().type(), dim-1);
for (size_t qp=0; qp<quad.size(); qp++) {
// Local position of the quadrature point
const Dune::FieldVector<double,dim-1>& quadPos = quad[qp].position();
const double integrationElement = nIt.intersectionGlobal().integrationElement(quadPos);
// Evaluate function
Dune::FieldVector<double,dim> localPressure(0);

Oliver Sander
committed
for (size_t i=0; i<baseSet.size(); i++)
localPressure.axpy(baseSet[i].evaluateFunction(0,quadPos),
pressure[dgIndexSet(*eIt,nIt.numberInSelf())+i]);

Oliver Sander
committed
// Sum up the total force
outputForce.axpy(quad[qp].weight()*integrationElement, localPressure);
// Sum up the total torque \int (x - x_0) \times f dx
Dune::FieldVector<double,dim> worldPos = nIt.intersectionGlobal().global(quadPos);
outputTorque.axpy(quad[qp].weight()*integrationElement,
crossProduct(worldPos - crossSection.r, localPressure));
}
}
}
outputForce -= resultantForce;
outputTorque -= resultantTorque;
assert( outputForce.two_norm() < 1e-6 );
assert( outputTorque.two_norm() < 1e-6 );
// std::cout << "Output force: " << outputForce << std::endl;
// std::cout << "Output torque: " << outputTorque << " " << resultantTorque[0]/outputTorque[0] << std::endl;

Oliver Sander
committed
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
template <class GridType>
void computeAverageInterface(const BoundaryPatch<GridType>& interface,
const Dune::BlockVector<Dune::FieldVector<double,GridType::dimension> > deformation,
Configuration& average)
{
using namespace Dune;
typedef typename GridType::template Codim<0>::LevelIterator ElementIterator;
typedef typename GridType::template Codim<0>::Entity EntityType;
typedef typename EntityType::LevelIntersectionIterator NeighborIterator;
const GridType& grid = interface.getGrid();
const int level = interface.level();
const typename GridType::Traits::LevelIndexSet& indexSet = grid.levelIndexSet(level);
const int dim = GridType::dimension;
// ///////////////////////////////////////////
// Initialize output configuration
// ///////////////////////////////////////////
average.r = 0;
double interfaceArea = 0;
FieldMatrix<double,dim,dim> deformationGradient(0);
// ///////////////////////////////////////////
// Loop and integrate over the interface
// ///////////////////////////////////////////
ElementIterator eIt = grid.template lbegin<0>(level);
ElementIterator eEndIt = grid.template lend<0>(level);
for (; eIt!=eEndIt; ++eIt) {
NeighborIterator nIt = eIt->ilevelbegin();
NeighborIterator nEndIt = eIt->ilevelend();
for (; nIt!=nEndIt; ++nIt) {
if (!interface.contains(*eIt, nIt))
continue;
const typename NeighborIterator::Geometry& segmentGeometry = nIt.intersectionGlobal();
const ReferenceElement<double,dim>& refElement = ReferenceElements<double, dim>::general(eIt->geometry().type());
int nDofs = refElement.size(nIt.numberInSelf(),1,dim);
// Get quadrature rule
const QuadratureRule<double, dim-1>& quad = QuadratureRules<double, dim-1>::rule(segmentGeometry.type(), dim-1);
// Get set of shape functions on this segment
const typename LagrangeShapeFunctionSetContainer<double,double,dim>::value_type& sfs
= LagrangeShapeFunctions<double,double,dim>::general(eIt->geometry().type(),1);
/* Loop over all integration points */
for (int ip=0; ip<quad.size(); ip++) {
// Local position of the quadrature point
//const FieldVector<double,dim-1>& quadPos = quad[ip].position();
const FieldVector<double,dim> quadPos = nIt.intersectionSelfLocal().global(quad[ip].position());
const double integrationElement = segmentGeometry.integrationElement(quad[ip].position());
// Evaluate base functions
FieldVector<double,dim> posAtQuadPoint(0);
for(int i=0; i<eIt->geometry().corners(); i++) {
int idx = indexSet.template subIndex<dim>(*eIt, i);
// Deformation at the quadrature point
posAtQuadPoint.axpy(sfs[i].evaluateFunction(0,quadPos), deformation[idx]);
}
const FieldMatrix<double,dim,dim>& inv = eIt->geometry().jacobianInverseTransposed(quadPos);
/* Compute the weight of the current integration point */
double weight = quad[ip].weight() * integrationElement;
/**********************************************/
/* compute gradients of the shape functions */
/**********************************************/
std::vector<FieldVector<double, dim> > shapeGrads(eIt->geometry().corners());
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
for (int dof=0; dof<eIt->geometry().corners(); dof++) {
for (int i=0; i<dim; i++)
shapeGrads[dof][i] = sfs[dof].evaluateDerivative(0, i, quadPos);
// multiply with jacobian inverse
FieldVector<double,dim> tmp(0);
inv.umv(shapeGrads[dof], tmp);
shapeGrads[dof] = tmp;
//std::cout << "Gradient " << dof << ": " << shape_grads[dof] << std::endl;
}
/****************************************************/
// The deformation gradient of the deformation
// in formulas: F(i,j) = $\partial \phi_i / \partial x_j$
// or F(i,j) = Id + $\partial u_i / \partial x_j$
/****************************************************/
FieldMatrix<double, dim, dim> F;
for (int i=0; i<dim; i++) {
for (int j=0; j<dim; j++) {
F[i][j] = (i==j) ? 1 : 0;
for (int k=0; k<eIt->geometry().corners(); k++)
F[i][j] += deformation[indexSet.template subIndex<dim>(*eIt, k)][i]*shapeGrads[k][j];
}
}
interfaceArea += quad[ip].weight() * integrationElement;
average.r.axpy(quad[ip].weight() * integrationElement, posAtQuadPoint);
F *= quad[ip].weight();
deformationGradient += F;
}
}
}
// average deformation of the interface is the integral over its
// deformation divided by its area
average.r /= interfaceArea;
// average deformation is the integral over the deformation gradient
// divided by its area
deformationGradient /= interfaceArea;
std::cout << "deformationGradient: " << std::endl << deformationGradient << std::endl;
// Get the rotational part of the deformation gradient by performing a
// polar composition.
FieldVector<double,dim> W;
FieldMatrix<double,dim,dim> VT;
// returns a decomposition U W VT, where U is returned in the first argument
svdcmp<double,dim,dim>(deformationGradient, W, VT);
deformationGradient.rightmultiply(VT);
// deformationGradient now contains the orthogonal part of the polar decomposition
assert( std::abs(1-deformationGradient.determinant()) < 1e-3);
}
#endif