From ad24249f24485ccb357e93e83003146393b0c5c8 Mon Sep 17 00:00:00 2001 From: Oliver Sander <sander@igpm.rwth-aachen.de> Date: Tue, 14 Dec 2010 16:46:11 +0000 Subject: [PATCH] Make computeAveragePressure take the centerOfTorque as a FieldVector. It was a RigidBodyMotion before, but the rotation part of that was ignored. Semantically it is just a position, hence a FieldVector is better. [[Imported from SVN: r6598]] --- dirneucoupling.cc | 4 ++-- dune/gfe/averageinterface.hh | 6 +++--- 2 files changed, 5 insertions(+), 5 deletions(-) diff --git a/dirneucoupling.cc b/dirneucoupling.cc index 31bffb14..55dfb41f 100644 --- a/dirneucoupling.cc +++ b/dirneucoupling.cc @@ -392,7 +392,7 @@ int main (int argc, char *argv[]) try // Using that index 0 is always the left boundary for a uniformly refined OneDGrid computeAveragePressure<GridType::LevelGridView>(resultantForce, resultantTorque, interfaceBoundary[grid.maxLevel()], - rodX[0], + rodX[0].r, neumannValues); rhs3d = 0; @@ -520,7 +520,7 @@ int main (int argc, char *argv[]) try // Using that index 0 is always the left boundary for a uniformly refined OneDGrid computeAveragePressure<GridType::LevelGridView>(resultantForce, resultantTorque, interfaceBoundary[grid.maxLevel()], - rodX[0], + rodX[0].r, neumannValues); rhs3d = 0; diff --git a/dune/gfe/averageinterface.hh b/dune/gfe/averageinterface.hh index 944633a3..ef479816 100644 --- a/dune/gfe/averageinterface.hh +++ b/dune/gfe/averageinterface.hh @@ -472,7 +472,7 @@ template <class GridView> void computeAveragePressure(const Dune::FieldVector<double,GridView::dimension>& resultantForce, const Dune::FieldVector<double,GridView::dimension>& resultantTorque, const BoundaryPatchBase<GridView>& interface, - const RigidBodyMotion<3>& crossSection, + const Dune::FieldVector<double,GridView::dimension>& centerOfTorque, Dune::BlockVector<Dune::FieldVector<double, GridView::dimension> >& pressure) { const GridView& gridView = interface.gridView(); @@ -567,7 +567,7 @@ void computeAveragePressure(const Dune::FieldVector<double,GridView::dimension>& phi_i[j] = shapeFunctionValues[i]; mu_tilde[i][j].axpy(quad[qp].weight() * integrationElement, - crossProduct(Dune::FieldVector<double,dim>(worldPos-crossSection.r), phi_i)); + crossProduct(Dune::FieldVector<double,dim>(worldPos-centerOfTorque), phi_i)); } @@ -653,7 +653,7 @@ void computeAveragePressure(const Dune::FieldVector<double,GridView::dimension>& Dune::FieldVector<double,3> outputForce(0), outputTorque(0); - computeTotalForceAndTorque(interface, pressure, crossSection.r, outputForce, outputTorque); + computeTotalForceAndTorque(interface, pressure, centerOfTorque, outputForce, outputTorque); outputForce -= resultantForce; outputTorque -= resultantTorque; -- GitLab