Skip to content
Snippets Groups Projects
Commit ceac2e3b authored by Oliver Sander's avatar Oliver Sander Committed by sander@PCPOOL.MI.FU-BERLIN.DE
Browse files

Adapt neumannassembler to use the dune-localfunctions shape functions.

This appears to be the last dependency on dune-disc of this module.
Hence I remove 'dune-disc' in dune.module.

[[Imported from SVN: r3041]]
parent 8ea5b799
No related branches found
No related tags found
No related merge requests found
#ifndef NEUMANN_ASSEMBLER_HH
#define NEUMANN_ASSEMBLER_HH
#if HAVE_DUNE_DISC
/** \file
\brief Contains a method which assembles the right hand side term for Neumann boundary data.
*/
......@@ -14,10 +12,8 @@
#include <dune/istl/matrix.hh>
#include <dune/disc/shapefunctions/lagrangeshapefunctions.hh>
#include <dune/ag-common/dgindexset.hh>
#include "boundarypatch.hh"
#include <dune/ag-common/boundarypatch.hh>
#include <dune/ag-common/functionspacebases/p1nodalbasis.hh>
template <class GridView, class VectorType>
void assembleAndAddNeumannTerm(const BoundaryPatchBase<GridView>& neumannBoundary,
......@@ -27,6 +23,9 @@ void assembleAndAddNeumannTerm(const BoundaryPatchBase<GridView>& neumannBoundar
const GridView& gridView = neumannBoundary.gridView();
const int dim = GridView::dimension;
// we need a P1 nodal basis as the shape function factory
P1NodalBasis<GridView,double> p1Basis(gridView);
typename BoundaryPatchBase<GridView>::iterator it = neumannBoundary.begin();
typename BoundaryPatchBase<GridView>::iterator endIt = neumannBoundary.end();
......@@ -47,9 +46,8 @@ void assembleAndAddNeumannTerm(const BoundaryPatchBase<GridView>& neumannBoundar
const int quadOrder = 4;
const Dune::QuadratureRule<double, dim-1>& quad = Dune::QuadratureRules<double, dim-1>::rule(segmentGeometry.type(), quadOrder);
// Get set of shape functions on this segment
const Dune::LagrangeShapeFunctionSet<double,double,dim-1>& sfs
= Dune::LagrangeShapeFunctions<double,double,dim-1>::general(segmentGeometry.type(),1);
// Get set of shape functions on this element
const typename P1NodalBasis<GridView,double>::LocalFiniteElement& localFiniteElement = p1Basis.getLocalFiniteElement(*it->inside());
/* Loop over all integration points */
for (int ip=0; ip<quad.size(); ip++) {
......@@ -59,12 +57,12 @@ void assembleAndAddNeumannTerm(const BoundaryPatchBase<GridView>& neumannBoundar
const double integrationElement = segmentGeometry.integrationElement(quadPos);
// Position of the quadrature point within the element
const Dune::FieldVector<double,dim> elementQuadPos = it->geometryInInside().global(quadPos);
// Evaluate base functions
double v[nDofs];
for (int i=0; i<nDofs; i++)
v[i] = sfs[i].evaluateFunction(0,quadPos);
std::vector<Dune::FieldVector<double,1> > v(nDofs);
localFiniteElement.localBasis().evaluateFunction(elementQuadPos,v);
for(int i=0; i<nDofs; i++)
for (int j=0; j<=i; j++ )
......@@ -98,8 +96,4 @@ void assembleAndAddNeumannTerm(const BoundaryPatchBase<GridView>& neumannBoundar
}
#else
#warning neumannassembler.hh was included but dune-disc is not present !
#endif // end of #if HAVE_DUNE_DISC
#endif
Module: ag-common
Depends: dune-common dune-grid dune-istl dune-localfunctions
Suggests: dune-disc
Suggests:
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment