displacementInterpolationMotionSolver.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) 2011-2012 OpenFOAM Foundation 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 \*---------------------------------------------------------------------------*/ 00025 00026 #include "displacementInterpolationMotionSolver.H" 00027 #include "addToRunTimeSelectionTable.H" 00028 #include "SortableList.H" 00029 #include "IOList.H" 00030 #include "Tuple2.H" 00031 #include "mapPolyMesh.H" 00032 #include "interpolateXY.H" 00033 00034 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // 00035 00036 namespace Foam 00037 { 00038 defineTypeNameAndDebug(displacementInterpolationMotionSolver, 0); 00039 00040 addToRunTimeSelectionTable 00041 ( 00042 motionSolver, 00043 displacementInterpolationMotionSolver, 00044 dictionary 00045 ); 00046 00047 template<> 00048 const word IOList<Tuple2<scalar, vector> >::typeName("scalarVectorTable"); 00049 } 00050 00051 00052 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // 00053 00054 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // 00055 00056 Foam::displacementInterpolationMotionSolver:: 00057 displacementInterpolationMotionSolver 00058 ( 00059 const polyMesh& mesh, 00060 const IOdictionary& dict 00061 ) 00062 : 00063 displacementMotionSolver(mesh, dict, typeName) 00064 { 00065 // Get zones and their interpolation tables for displacement 00066 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 00067 00068 List<Pair<word> > faceZoneToTable 00069 ( 00070 coeffDict().lookup("interpolationTables") 00071 ); 00072 00073 const faceZoneMesh& fZones = mesh.faceZones(); 00074 00075 times_.setSize(fZones.size()); 00076 displacements_.setSize(fZones.size()); 00077 00078 forAll(faceZoneToTable, i) 00079 { 00080 const word& zoneName = faceZoneToTable[i][0]; 00081 label zoneI = fZones.findZoneID(zoneName); 00082 00083 if (zoneI == -1) 00084 { 00085 FatalErrorIn 00086 ( 00087 "displacementInterpolationMotionSolver::" 00088 "displacementInterpolationMotionSolver(const polyMesh&," 00089 "Istream&)" 00090 ) << "Cannot find zone " << zoneName << endl 00091 << "Valid zones are " << mesh.faceZones().names() 00092 << exit(FatalError); 00093 } 00094 00095 const word& tableName = faceZoneToTable[i][1]; 00096 00097 IOList<Tuple2<scalar, vector> > table 00098 ( 00099 IOobject 00100 ( 00101 tableName, 00102 mesh.time().constant(), 00103 "tables", 00104 mesh, 00105 IOobject::MUST_READ, 00106 IOobject::NO_WRITE, 00107 false 00108 ) 00109 ); 00110 00111 // Copy table 00112 times_[zoneI].setSize(table.size()); 00113 displacements_[zoneI].setSize(table.size()); 00114 00115 forAll(table, j) 00116 { 00117 times_[zoneI][j] = table[j].first(); 00118 displacements_[zoneI][j] = table[j].second(); 00119 } 00120 } 00121 00122 00123 00124 // Sort points into bins according to position relative to faceZones 00125 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 00126 // Done in all three directions. 00127 00128 for (direction dir = 0; dir < vector::nComponents; dir++) 00129 { 00130 // min and max coordinates of all faceZones 00131 SortableList<scalar> zoneCoordinates(2*faceZoneToTable.size()); 00132 00133 forAll(faceZoneToTable, i) 00134 { 00135 const word& zoneName = faceZoneToTable[i][0]; 00136 const faceZone& fz = fZones[zoneName]; 00137 00138 scalar minCoord = VGREAT; 00139 scalar maxCoord = -VGREAT; 00140 00141 forAll(fz().meshPoints(), localI) 00142 { 00143 label pointI = fz().meshPoints()[localI]; 00144 const scalar coord = points0()[pointI][dir]; 00145 minCoord = min(minCoord, coord); 00146 maxCoord = max(maxCoord, coord); 00147 } 00148 00149 zoneCoordinates[2*i] = returnReduce(minCoord, minOp<scalar>()); 00150 zoneCoordinates[2*i+1] = returnReduce(maxCoord, maxOp<scalar>()); 00151 00152 if (debug) 00153 { 00154 Pout<< "direction " << dir << " : " 00155 << "zone " << zoneName 00156 << " ranges from coordinate " << zoneCoordinates[2*i] 00157 << " to " << zoneCoordinates[2*i+1] 00158 << endl; 00159 } 00160 } 00161 zoneCoordinates.sort(); 00162 00163 // Slightly tweak min and max face zone so points sort within 00164 zoneCoordinates[0] -= SMALL; 00165 zoneCoordinates.last() += SMALL; 00166 00167 // Check if we have static min and max mesh bounds 00168 const scalarField meshCoords(points0().component(dir)); 00169 00170 scalar minCoord = gMin(meshCoords); 00171 scalar maxCoord = gMax(meshCoords); 00172 00173 if (debug) 00174 { 00175 Pout<< "direction " << dir << " : " 00176 << "mesh ranges from coordinate " << minCoord << " to " 00177 << maxCoord << endl; 00178 } 00179 00180 // Make copy of zoneCoordinates; include min and max of mesh 00181 // if necessary. Mark min and max with zoneI=-1. 00182 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 00183 00184 labelList& rangeZone = rangeToZone_[dir]; 00185 labelListList& rangePoints = rangeToPoints_[dir]; 00186 List<scalarField>& rangeWeights = rangeToWeights_[dir]; 00187 00188 scalarField rangeToCoord(zoneCoordinates.size()); 00189 rangeZone.setSize(zoneCoordinates.size()); 00190 label rangeI = 0; 00191 00192 if (minCoord < zoneCoordinates[0]) 00193 { 00194 label sz = rangeZone.size(); 00195 rangeToCoord.setSize(sz+1); 00196 rangeZone.setSize(sz+1); 00197 rangeToCoord[rangeI] = minCoord-SMALL; 00198 rangeZone[rangeI] = -1; 00199 00200 if (debug) 00201 { 00202 Pout<< "direction " << dir << " : " 00203 << "range " << rangeI << " at coordinate " 00204 << rangeToCoord[rangeI] << " from min of mesh " 00205 << rangeZone[rangeI] << endl; 00206 } 00207 rangeI = 1; 00208 } 00209 forAll(zoneCoordinates, i) 00210 { 00211 rangeToCoord[rangeI] = zoneCoordinates[i]; 00212 rangeZone[rangeI] = zoneCoordinates.indices()[i]/2; 00213 00214 if (debug) 00215 { 00216 Pout<< "direction " << dir << " : " 00217 << "range " << rangeI << " at coordinate " 00218 << rangeToCoord[rangeI] 00219 << " from zone " << rangeZone[rangeI] << endl; 00220 } 00221 rangeI++; 00222 } 00223 if (maxCoord > zoneCoordinates.last()) 00224 { 00225 label sz = rangeToCoord.size(); 00226 rangeToCoord.setSize(sz+1); 00227 rangeZone.setSize(sz+1); 00228 rangeToCoord[sz] = maxCoord+SMALL; 00229 rangeZone[sz] = -1; 00230 00231 if (debug) 00232 { 00233 Pout<< "direction " << dir << " : " 00234 << "range " << rangeI << " at coordinate " 00235 << rangeToCoord[sz] << " from max of mesh " 00236 << rangeZone[sz] << endl; 00237 } 00238 } 00239 00240 00241 // Sort the points 00242 // ~~~~~~~~~~~~~~~ 00243 00244 // Count all the points inbetween rangeI and rangeI+1 00245 labelList nRangePoints(rangeToCoord.size(), 0); 00246 00247 forAll(meshCoords, pointI) 00248 { 00249 label rangeI = findLower(rangeToCoord, meshCoords[pointI]); 00250 00251 if (rangeI == -1 || rangeI == rangeToCoord.size()-1) 00252 { 00253 FatalErrorIn 00254 ( 00255 "displacementInterpolationMotionSolver::" 00256 "displacementInterpolationMotionSolver" 00257 "(const polyMesh&, Istream&)" 00258 ) << "Did not find point " << points0()[pointI] 00259 << " coordinate " << meshCoords[pointI] 00260 << " in ranges " << rangeToCoord 00261 << abort(FatalError); 00262 } 00263 nRangePoints[rangeI]++; 00264 } 00265 00266 if (debug) 00267 { 00268 for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++) 00269 { 00270 // Get the two zones bounding the range 00271 Pout<< "direction " << dir << " : " 00272 << "range from " << rangeToCoord[rangeI] 00273 << " to " << rangeToCoord[rangeI+1] 00274 << " contains " << nRangePoints[rangeI] 00275 << " points." << endl; 00276 } 00277 } 00278 00279 // Sort 00280 rangePoints.setSize(nRangePoints.size()); 00281 rangeWeights.setSize(nRangePoints.size()); 00282 forAll(rangePoints, rangeI) 00283 { 00284 rangePoints[rangeI].setSize(nRangePoints[rangeI]); 00285 rangeWeights[rangeI].setSize(nRangePoints[rangeI]); 00286 } 00287 nRangePoints = 0; 00288 forAll(meshCoords, pointI) 00289 { 00290 label rangeI = findLower(rangeToCoord, meshCoords[pointI]); 00291 label& nPoints = nRangePoints[rangeI]; 00292 rangePoints[rangeI][nPoints] = pointI; 00293 rangeWeights[rangeI][nPoints] = 00294 (meshCoords[pointI]-rangeToCoord[rangeI]) 00295 / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]); 00296 nPoints++; 00297 } 00298 } 00299 } 00300 00301 00302 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // 00303 00304 Foam::displacementInterpolationMotionSolver:: 00305 ~displacementInterpolationMotionSolver() 00306 {} 00307 00308 00309 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // 00310 00311 Foam::tmp<Foam::pointField> 00312 Foam::displacementInterpolationMotionSolver::curPoints() const 00313 { 00314 if (mesh().nPoints() != points0().size()) 00315 { 00316 FatalErrorIn 00317 ( 00318 "displacementInterpolationMotionSolver::curPoints() const" 00319 ) << "The number of points in the mesh seems to have changed." << endl 00320 << "In constant/polyMesh there are " << points0().size() 00321 << " points; in the current mesh there are " << mesh().nPoints() 00322 << " points." << exit(FatalError); 00323 } 00324 00325 tmp<pointField> tcurPoints(new pointField(points0())); 00326 pointField& curPoints = tcurPoints(); 00327 00328 // Interpolate the diplacement of the face zones. 00329 vectorField zoneDisp(displacements_.size(), vector::zero); 00330 forAll(zoneDisp, zoneI) 00331 { 00332 if (times_[zoneI].size()) 00333 { 00334 zoneDisp[zoneI] = interpolateXY 00335 ( 00336 mesh().time().value(), 00337 times_[zoneI], 00338 displacements_[zoneI] 00339 ); 00340 } 00341 } 00342 if (debug) 00343 { 00344 Pout<< "Zone displacements:" << zoneDisp << endl; 00345 } 00346 00347 00348 // Interpolate the point location 00349 for (direction dir = 0; dir < vector::nComponents; dir++) 00350 { 00351 const labelList& rangeZone = rangeToZone_[dir]; 00352 const labelListList& rangePoints = rangeToPoints_[dir]; 00353 const List<scalarField>& rangeWeights = rangeToWeights_[dir]; 00354 00355 for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++) 00356 { 00357 const labelList& rPoints = rangePoints[rangeI]; 00358 const scalarField& rWeights = rangeWeights[rangeI]; 00359 00360 // Get the two zones bounding the range 00361 label minZoneI = rangeZone[rangeI]; 00362 //vector minDisp = 00363 // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]); 00364 scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]); 00365 label maxZoneI = rangeZone[rangeI+1]; 00366 //vector maxDisp = 00367 // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]); 00368 scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]); 00369 00370 forAll(rPoints, i) 00371 { 00372 label pointI = rPoints[i]; 00373 scalar w = rWeights[i]; 00374 //curPoints[pointI] += (1.0-w)*minDisp+w*maxDisp; 00375 curPoints[pointI][dir] += (1.0-w)*minDisp+w*maxDisp; 00376 } 00377 } 00378 } 00379 return tcurPoints; 00380 } 00381 00382 00383 // ************************************************************************* //
