FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

edgeVertex.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Description
00025 
00026 \*---------------------------------------------------------------------------*/
00027 
00028 #include "edgeVertex.H"
00029 #include <meshTools/meshTools.H>
00030 #include <dynamicMesh/refineCell.H>
00031 
00032 
00033 // * * * * * * * * * * * * * * * Static Functions  * * * * * * * * * * * * * //
00034 
00035 
00036 // Update stored refine list using map
00037 void Foam::edgeVertex::updateLabels
00038 (
00039     const labelList& map,
00040     List<refineCell>& refCells
00041 )
00042 {
00043     label newRefI = 0;
00044 
00045     forAll(refCells, refI)
00046     {
00047         const refineCell& refCell = refCells[refI];
00048 
00049         label oldCellI = refCell.cellNo();
00050 
00051         label newCellI = map[oldCellI];
00052 
00053         if (newCellI != -1)
00054         {
00055             refCells[newRefI++] = refineCell(newCellI, refCell.direction());
00056         }
00057     }
00058     refCells.setSize(newRefI);
00059 }
00060 
00061 
00062 // Update stored cell numbers using map.
00063 // Do in two passes to prevent allocation if nothing changed.
00064 void Foam::edgeVertex::updateLabels
00065 (
00066     const labelList& map,
00067     Map<label>& cellPairs
00068 )
00069 {
00070     // Iterate over map to see if anything changed
00071     bool changed = false;
00072 
00073     for
00074     (
00075         Map<label>::const_iterator iter = cellPairs.begin();
00076         iter != cellPairs.end();
00077         ++iter
00078     )
00079     {
00080         label newMaster = map[iter.key()];
00081 
00082         label newSlave = -1;
00083 
00084         if (iter() != -1)
00085         {
00086             newSlave = map[iter()];
00087         }
00088 
00089         if ((newMaster != iter.key()) || (newSlave != iter()))
00090         {
00091             changed = true;
00092 
00093             break;
00094         }
00095     }
00096 
00097     // Relabel (use second Map to prevent overlapping)
00098     if (changed)
00099     {
00100         Map<label> newCellPairs(2*cellPairs.size());
00101 
00102         for
00103         (
00104             Map<label>::const_iterator iter = cellPairs.begin();
00105             iter != cellPairs.end();
00106             ++iter
00107         )
00108         {
00109             label newMaster = map[iter.key()];
00110 
00111             label newSlave = -1;
00112 
00113             if (iter() != -1)
00114             {
00115                 newSlave = map[iter()];
00116             }
00117 
00118             if (newMaster == -1)
00119             {
00120                 WarningIn
00121                 (
00122                     "edgeVertex::updateLabels(const labelList&, "
00123                     "Map<label>&)"
00124                 )   << "master cell:" << iter.key()
00125                     << " has disappeared" << endl;
00126             }
00127             else
00128             {
00129                 newCellPairs.insert(newMaster, newSlave);
00130             }
00131         }
00132 
00133         cellPairs = newCellPairs;
00134     }
00135 }
00136 
00137 
00138 // Update stored cell numbers using map.
00139 // Do in two passes to prevent allocation if nothing changed.
00140 void Foam::edgeVertex::updateLabels
00141 (
00142     const labelList& map,
00143     labelHashSet& cells
00144 )
00145 {
00146     // Iterate over map to see if anything changed
00147     bool changed = false;
00148 
00149     for
00150     (
00151         labelHashSet::const_iterator iter = cells.begin();
00152         iter != cells.end();
00153         ++iter
00154     )
00155     {
00156         label newCellI = map[iter.key()];
00157 
00158         if (newCellI != iter.key())
00159         {
00160             changed = true;
00161 
00162             break;
00163         }
00164     }
00165 
00166     // Relabel (use second Map to prevent overlapping)
00167     if (changed)
00168     {
00169         labelHashSet newCells(2*cells.size());
00170 
00171         for
00172         (
00173             labelHashSet::const_iterator iter = cells.begin();
00174             iter != cells.end();
00175             ++iter
00176         )
00177         {
00178             label newCellI = map[iter.key()];
00179 
00180             if (newCellI != -1)
00181             {
00182                 newCells.insert(newCellI);
00183             }
00184         }
00185 
00186         cells = newCells;
00187     }
00188 }
00189 
00190 
00191 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00192 
00193 Foam::point Foam::edgeVertex::coord
00194 (
00195     const primitiveMesh& mesh,
00196     const label cut,
00197     const scalar weight
00198 )
00199 {
00200     const pointField& pts = mesh.points();
00201 
00202     if (isEdge(mesh, cut))
00203     {
00204         const edge& e = mesh.edges()[getEdge(mesh, cut)];
00205 
00206         return weight*pts[e.end()] + (1-weight)*pts[e.start()];
00207     }
00208     else
00209     {
00210         return pts[getVertex(mesh, cut)];
00211     }
00212 }
00213 
00214 
00215 Foam::label Foam::edgeVertex::cutPairToEdge
00216 (
00217     const primitiveMesh& mesh,
00218     const label cut0,
00219     const label cut1
00220 )
00221 {
00222     if (!isEdge(mesh, cut0) && !isEdge(mesh, cut1))
00223     {
00224         return meshTools::findEdge
00225         (
00226             mesh,
00227             getVertex(mesh, cut0),
00228             getVertex(mesh, cut1)
00229         );
00230     }
00231     else
00232     {
00233         return -1;
00234     }
00235 }
00236 
00237 
00238 Foam::Ostream& Foam::edgeVertex::writeCut
00239 (
00240     Ostream& os,
00241     const label cut,
00242     const scalar weight
00243 ) const
00244 {
00245     if (isEdge(cut))
00246     {
00247         label edgeI = getEdge(cut);
00248 
00249         const edge& e = mesh().edges()[edgeI];
00250 
00251         os  << "edge:" << edgeI << e << ' ' << coord(cut, weight);
00252     }
00253     else
00254     {
00255         label vertI = getVertex(cut);
00256 
00257         os << "vertex:" << vertI << ' ' << coord(cut, weight);
00258     }
00259     return os;
00260 }
00261 
00262 
00263 Foam::Ostream& Foam::edgeVertex::writeCuts
00264 (
00265     Ostream& os,
00266     const labelList& cuts,
00267     const scalarField& weights
00268 ) const
00269 {
00270     forAll(cuts, cutI)
00271     {
00272         if (cutI > 0)
00273         {
00274             os << ' ';
00275         }
00276         writeCut(os, cuts[cutI], weights[cutI]);
00277     }
00278     return os;
00279 }
00280 
00281 
00282 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines