Skip to content
Snippets Groups Projects
Commit 177176c5 authored by mattijs's avatar mattijs Committed by Kutalmış Berçin
Browse files

BUG: particle: use correct celli. Fixes #1992

Was checking the old celli instead of the result of
re-finding the position. See also Foundation commit 50a965f8866683a81d79cbc7811af7333baf9d10.
parent 53feaf25
Branches
No related tags found
No related merge requests found
...@@ -6,7 +6,7 @@ ...@@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011-2017, 2020 OpenFOAM Foundation Copyright (C) 2011-2017, 2020 OpenFOAM Foundation
Copyright (C) 2018-2020 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
...@@ -415,21 +415,21 @@ void Foam::particle::locate ...@@ -415,21 +415,21 @@ void Foam::particle::locate
{ {
if (debug) if (debug)
{ {
Info << "Particle " << origId() << nl << FUNCTION_NAME << nl << endl; Pout<< "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
} }
celli_ = celli; celli_ = celli;
// Find the cell, if it has not been given // Find the cell, if it has not been given
if (celli < 0) if (celli_ < 0)
{ {
celli_ = mesh_.cellTree().findInside(position); celli_ = mesh_.cellTree().findInside(position);
} if (celli_ < 0)
if (celli < 0) {
{ FatalErrorInFunction
FatalErrorInFunction << "Cell not found for particle position " << position << "."
<< "Cell not found for particle position " << position << "." << exit(FatalError);
<< exit(FatalError); }
} }
const vector displacement = position - mesh_.cellCentres()[celli_]; const vector displacement = position - mesh_.cellCentres()[celli_];
...@@ -811,10 +811,10 @@ Foam::scalar Foam::particle::trackToStationaryTri ...@@ -811,10 +811,10 @@ Foam::scalar Foam::particle::trackToStationaryTri
if (debug) if (debug)
{ {
Info << "Step Fraction : " << stepFraction_ << endl; Pout << "Step Fraction : " << stepFraction_ << endl;
Info << "fraction*muH*detA : " << fraction*muH*detA << endl; Pout << "fraction*muH*detA : " << fraction*muH*detA << endl;
Info << "static muH : " << muH << endl; Pout << "static muH : " << muH << endl;
Info << "origId() : " << origId() << endl; Pout << "origId() : " << origId() << endl;
} }
// Accumulate displacement behind // Accumulate displacement behind
...@@ -864,7 +864,7 @@ Foam::scalar Foam::particle::trackToMovingTri ...@@ -864,7 +864,7 @@ Foam::scalar Foam::particle::trackToMovingTri
{ {
Pair<vector> o, b, v1, v2; Pair<vector> o, b, v1, v2;
movingTetGeometry(fraction, o, b, v1, v2); movingTetGeometry(fraction, o, b, v1, v2);
Info<< "Tet points o=" << o[0] << ", b=" << b[0] Pout<< "Tet points o=" << o[0] << ", b=" << b[0]
<< ", v1=" << v1[0] << ", v2=" << v2[0] << endl << ", v1=" << v1[0] << ", v2=" << v2[0] << endl
<< "Tet determinant = " << detA[0] << endl << "Tet determinant = " << detA[0] << endl
<< "Start local coordinates = " << y0[0] << endl; << "Start local coordinates = " << y0[0] << endl;
...@@ -895,10 +895,10 @@ Foam::scalar Foam::particle::trackToMovingTri ...@@ -895,10 +895,10 @@ Foam::scalar Foam::particle::trackToMovingTri
{ {
for (label i = 0; i < 4; ++ i) for (label i = 0; i < 4; ++ i)
{ {
Info<< (i ? " " : "Hit equation ") << i << " = " Pout<< (i ? " " : "Hit equation ") << i << " = "
<< hitEqn[i] << endl; << hitEqn[i] << endl;
} }
Info<< " DetA equation = " << detA << endl; Pout<< " DetA equation = " << detA << endl;
} }
// Calculate the hit fraction // Calculate the hit fraction
...@@ -927,12 +927,12 @@ Foam::scalar Foam::particle::trackToMovingTri ...@@ -927,12 +927,12 @@ Foam::scalar Foam::particle::trackToMovingTri
); );
const scalar detAH = detAEqn.value(mu[j]); const scalar detAH = detAEqn.value(mu[j]);
Info<< "Hit on tet face " << i << " at local coordinate " Pout<< "Hit on tet face " << i << " at local coordinate "
<< (std::isnormal(detAH) ? name(yH/detAH) : "???") << (std::isnormal(detAH) ? name(yH/detAH) : "???")
<< ", " << mu[j]*detA[0]*100 << "% of the " << ", " << mu[j]*detA[0]*100 << "% of the "
<< "way along the track" << endl; << "way along the track" << endl;
Info<< "derivative : " << hitEqn[i].derivative(mu[j]) << nl Pout<< "derivative : " << hitEqn[i].derivative(mu[j]) << nl
<< " coord " << j << " mu[j]: " << mu[j] << nl << " coord " << j << " mu[j]: " << mu[j] << nl
<< " hitEq " << i << endl; << " hitEq " << i << endl;
} }
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment