/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright held by original author \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM; if not, write to the Free Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA Description Create coupled match faces and add them to the cells \*---------------------------------------------------------------------------*/ #include "starMesh.H" #include "boolList.H" #include "pointHit.H" #include "IOmanip.H" #include "boundBox.H" #include "Map.H" #include "mathematicalConstants.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // void starMesh::createCoupleMatches() { // Loop through all couples and create intersection faces. Add all points // of intersection faces to the couple points lists. The numbering of // the list is set such that the list can be appended to the // existing points list // Estimate the number of cells affected by couple matches const label cellMapSize = min ( cellShapes_.size()/10, couples_.size()*2 ); // Store newly created faces for each cell Map > cellAddedFaces(cellMapSize); Map > cellRemovedFaces(cellMapSize); // In order to remove often allocation, remember the number of live points. // If you run out of space in point creation, increase it by the number of // couples (good scale) and resize at the end; label nLivePoints = points_.size(); const label infoJump = max(1000, couples_.size()/20); forAll (couples_, coupleI) { if (coupleI % infoJump == 0) { Info << "Doing couple " << coupleI << ". STAR couple ID: " << couples_[coupleI].coupleID() << endl; } // Initialise cell edges for master and slave cells const coupledFacePair& fp = couples_[coupleI]; const face& masterFace = cellFaces_[fp.masterCell()][fp.masterFace()]; const face& slaveFace = cellFaces_[fp.slaveCell()][fp.slaveFace()]; # ifdef DEBUG_COUPLE Info<< "coupleI: " << coupleI << endl << "masterFace: " << masterFace << endl << "master points: " << masterFace.points(points_) << endl << "slaveFace: " << slaveFace << endl << "slave points: " << slaveFace.points(points_) << endl << endl; # endif // check the angle of face area vectors scalar faceAreaAngle = mag ( -(masterFace.normal(points_) & slaveFace.normal(points_))/ (masterFace.mag(points_)*slaveFace.mag(points_) + VSMALL) ); if (faceAreaAngle < 0.94) { Info<< "Couple direction mismatch in the couple match " << coupleI << ". STAR couple ID: " << couples_[coupleI].coupleID() << endl << "The angle between face normals is " << Foam::acos(faceAreaAngle)/mathematicalConstant::pi*180 << " deg." << endl << "master cell: " << fp.masterCell() << " STAR number: " << starCellID_[fp.masterCell()] << " type: " << cellShapes_[fp.masterCell()].model().name() << " face: " << fp.masterFace() << endl << "slave cell : " << fp.slaveCell() << " STAR number: " << starCellID_[fp.slaveCell()] << " type: " << cellShapes_[fp.slaveCell()].model().name() << " face: " << fp.slaveFace() << endl; } // Deal with integral patches if (fp.integralMatch()) { // Master face is replaced by a set of slave faces Map >::iterator crfIter = cellRemovedFaces.find(fp.masterCell()); if (crfIter == cellRemovedFaces.end()) { cellRemovedFaces.insert ( fp.masterCell(), fp.masterFace() ); } else { crfIter().append(fp.masterFace()); } Map >::iterator cafIter = cellAddedFaces.find(fp.masterCell()); if (cafIter == cellAddedFaces.end()) { cellAddedFaces.insert ( fp.masterCell(), SLList(slaveFace.reverseFace()) ); } else { cafIter().append(slaveFace.reverseFace()); } } else { // Create cut faces, which replace both master and slave faces // Store newly created points SLList coupleFacePoints; // Master data edgeList masterEdges = masterFace.edges(); List > masterEdgePoints(masterEdges.size()); // Slave data edgeList slaveEdges = slaveFace.edges(); List > slaveEdgePoints(slaveEdges.size()); // Find common plane vector n = masterFace.normal(points_); n /= mag(n) + VSMALL; // Loop through all edges of the master face. For every edge, // intersect it with all edges of the cutting face. forAll (masterEdges, masterEdgeI) { const edge& curMasterEdge = masterEdges[masterEdgeI]; point P = points_[curMasterEdge.start()]; // get d and return it into plane vector d = curMasterEdge.vec(points_); d -= n*(n & d); # ifdef DEBUG_COUPLE_INTERSECTION Info << "curMasterEdge: " << curMasterEdge << endl << "P: " << P << endl << "d: " << d << endl; # endif // go through all slave edges and try to get an intersection. // The point is created along the original master edge rather // than its corrected direction. forAll (slaveEdges, slaveEdgeI) { const edge& curSlaveEdge = slaveEdges[slaveEdgeI]; point S = points_[curSlaveEdge.start()]; // get e and return it into plane vector e = curSlaveEdge.vec(points_); e -= n*(n & e); scalar det = -(e & (n ^ d)); # ifdef DEBUG_COUPLE_INTERSECTION Info << "curSlaveEdge: " << curSlaveEdge << endl << "S: " << S << endl << "e: " << e << endl; # endif if (mag(det) > SMALL) { // non-singular matrix. Look for intersection scalar beta = ((S - P) & (n ^ d))/det; # ifdef DEBUG_COUPLE_INTERSECTION Info << " beta: " << beta << endl; # endif if (beta > -smallMergeTol_ && beta < 1 + smallMergeTol_) { // slave intersection OK. Try master intersection scalar alpha = (((S - P) & d) + beta*(d & e))/magSqr(d); # ifdef DEBUG_COUPLE_INTERSECTION Info << " alpha: " << alpha << endl; # endif if ( alpha > -smallMergeTol_ && alpha < 1 + smallMergeTol_ ) { // intersection of non-parallel edges # ifdef DEBUG_COUPLE_INTERSECTION Info<< "intersection of non-parallel edges" << endl; # endif // check for insertion of start-end // points in the middle of the other // edge if (alpha < smallMergeTol_) { // inserting the start of master edge if ( beta > smallMergeTol_ && beta < 1 - smallMergeTol_ ) { slaveEdgePoints[slaveEdgeI].append ( curMasterEdge.start() ); } } else if (alpha > 1 - smallMergeTol_) { // inserting the end of master edge if ( beta > smallMergeTol_ && beta < 1 - smallMergeTol_ ) { slaveEdgePoints[slaveEdgeI].append ( curMasterEdge.end() ); } } else if (beta < smallMergeTol_) { // inserting the start of the slave edge if ( alpha > smallMergeTol_ && alpha < 1 - smallMergeTol_ ) { masterEdgePoints[masterEdgeI].append ( curSlaveEdge.start() ); } } else if (beta > 1 - smallMergeTol_) { // inserting the start of the slave edge if ( alpha > smallMergeTol_ && alpha < 1 - smallMergeTol_ ) { masterEdgePoints[masterEdgeI].append ( curSlaveEdge.end() ); } } else { masterEdgePoints[masterEdgeI].append ( nLivePoints + coupleFacePoints.size() ); slaveEdgePoints[slaveEdgeI].append ( nLivePoints + coupleFacePoints.size() ); # ifdef DEBUG_COUPLE_INTERSECTION Info<< "regular intersection. " << "Adding point: " << coupleFacePoints.size() << " which is " << P + alpha*curMasterEdge.vec(points_) << endl; # endif // A new point is created. Warning: // using original edge for accuracy. // coupleFacePoints.append (P + alpha*curMasterEdge.vec(points_)); } } } } else { // Add special cases, for intersection of two // parallel line Warning. Here, typically, no new // points will be created. Either one or two of // the slave edge points need to be added to the // master edge and vice versa. The problem is that // no symmetry exists, i.e. both operations needs // to be done separately for both master and slave // side. // Master side // check if the first or second point of slave edge is // on the master edge vector ps = S - P; bool colinear = false; if (mag(ps) < SMALL) { // colinear because P and S are the same point colinear = true; } else if ( (ps & d)/(mag(ps)*mag(d)) > 1.0 - smallMergeTol_ ) { // colinear because ps and d are parallel colinear = true; } if (colinear) { scalar alpha1 = (ps & d)/magSqr(d); if ( alpha1 > -smallMergeTol_ && alpha1 < 1 + smallMergeTol_ ) { # ifdef DEBUG_COUPLE_INTERSECTION Info<< "adding irregular master " << "intersection1: " << points_[slaveEdges[slaveEdgeI].start()] << endl; # endif masterEdgePoints[masterEdgeI].append ( slaveEdges[slaveEdgeI].start() ); } scalar alpha2 = ((ps + e) & d)/magSqr(d); if ( alpha2 > -smallMergeTol_ && alpha2 < 1 + smallMergeTol_ ) { # ifdef DEBUG_COUPLE_INTERSECTION Info<< "adding irregular master " << "intersection2: " << points_[slaveEdges[slaveEdgeI].end()] << endl; # endif masterEdgePoints[masterEdgeI].append ( slaveEdges[slaveEdgeI].end() ); } // Slave side // check if the first or second point of // master edge is on the slave edge vector sp = P - S; scalar beta1 = (sp & e)/magSqr(e); # ifdef DEBUG_COUPLE_INTERSECTION Info << "P: " << P << " S: " << S << " d: " << d << " e: " << e << " sp: " << sp << " beta1: " << beta1 << endl; # endif if ( beta1 > -smallMergeTol_ && beta1 < 1 + smallMergeTol_ ) { # ifdef DEBUG_COUPLE_INTERSECTION Info<< "adding irregular slave " << "intersection1: " << points_[masterEdges[masterEdgeI].start()] << endl; # endif slaveEdgePoints[slaveEdgeI].append ( masterEdges[masterEdgeI].start() ); } scalar beta2 = ((sp + d) & e)/magSqr(e); if ( beta2 > -smallMergeTol_ && beta2 < 1 + smallMergeTol_ ) { # ifdef DEBUG_COUPLE_INTERSECTION Info << "adding irregular slave " << "intersection2: " << points_[masterEdges[masterEdgeI].end()] << endl; # endif slaveEdgePoints[slaveEdgeI].append ( masterEdges[masterEdgeI].end() ); } } // end of colinear } // end of singular intersection } // end of slave edges } // end of master edges # ifdef DEBUG_COUPLE_INTERSECTION Info << "additional slave edge points: " << endl; forAll (slaveEdgePoints, edgeI) { Info << "edge: " << edgeI << ": " << slaveEdgePoints[edgeI] << endl; } # endif // Add new points if (nLivePoints + coupleFacePoints.size() >= points_.size()) { // increase the size of the points list Info << "Resizing points list" << endl; points_.setSize(points_.size() + couples_.size()); } for ( SLList::iterator coupleFacePointsIter = coupleFacePoints.begin(); coupleFacePointsIter != coupleFacePoints.end(); ++coupleFacePointsIter ) { points_[nLivePoints] = coupleFacePointsIter(); nLivePoints++; } // edge intersection finished // Creating new master side // count the number of additional points for face label nAdditionalMasterPoints = 0; forAll (masterEdgePoints, edgeI) { nAdditionalMasterPoints += masterEdgePoints[edgeI].size(); } face tmpMasterFace ( masterFace.size() + nAdditionalMasterPoints ); label nTmpMasterLabels = 0; # ifdef DEBUG_COUPLE_INTERSECTION Info << "masterFace: " << masterFace << endl << "nAdditionalMasterPoints: " << nAdditionalMasterPoints << endl; # endif forAll (masterEdges, masterEdgeI) { // Insert the starting point of the edge tmpMasterFace[nTmpMasterLabels] = masterEdges[masterEdgeI].start(); nTmpMasterLabels++; // get reference to added points of current edge const SLList