Test-spline.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "argList.H"
27 
28 #include "vector.H"
29 #include "IFstream.H"
30 
31 #include "BSpline.H"
32 #include "CatmullRomSpline.H"
33 
34 using namespace Foam;
35 
36 inline Ostream& printPoint(Ostream& os, const point& p)
37 {
38  os << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
39  return os;
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 // Main program:
45 
46 int main(int argc, char *argv[])
47 {
49  argList::validArgs.insert("file .. fileN");
50  argList::addBoolOption("B", "B-Spline implementation");
51  argList::addBoolOption("CMR", "catmull-rom spline (default)");
53  (
54  "n",
55  "INT",
56  "number of segments for evaluation - default 20"
57  );
58 
59  argList args(argc, argv, false, true);
60 
61  if (args.size() <= 1)
62  {
63  args.printUsage();
64  }
65 
66  bool useBSpline = args.optionFound("B");
67  bool useCatmullRom = args.optionFound("CMR");
68  label nSeg = args.optionLookupOrDefault<label>("n", 20);
69 
70  if (!useCatmullRom && !useBSpline)
71  {
72  Info<<"defaulting to Catmull-Rom spline" << endl;
73  useCatmullRom = true;
74  }
75 
76  for (label argI=1; argI < args.size(); ++argI)
77  {
78  const string& srcFile = args[argI];
79  Info<< nl << "reading " << srcFile << nl;
80  IFstream ifs(srcFile);
81 
82  List<pointField> pointFields(ifs);
83 
84 
85  forAll(pointFields, splineI)
86  {
87  Info<<"\n# points:" << endl;
88  forAll(pointFields[splineI], ptI)
89  {
90  printPoint(Info, pointFields[splineI][ptI]);
91  }
92 
93  if (useBSpline)
94  {
95  BSpline spl(pointFields[splineI]);
96 
97  Info<< nl << "# B-Spline" << endl;
98 
99  for (label segI = 0; segI <= nSeg; ++segI)
100  {
101  scalar lambda = scalar(segI)/scalar(nSeg);
103  }
104  }
105 
106  if (useCatmullRom)
107  {
108  CatmullRomSpline spl(pointFields[splineI]);
109 
110  Info<< nl <<"# Catmull-Rom" << endl;
111 
112  for (label segI = 0; segI <= nSeg; ++segI)
113  {
114  scalar lambda = scalar(segI)/scalar(nSeg);
116  }
117  }
118 
119  {
120  polyLine pl(pointFields[splineI]);
121 
122  Info<< nl <<"# polyList" << endl;
123 
124  for (label segI = 0; segI <= nSeg; ++segI)
125  {
126  scalar lambda = scalar(segI)/scalar(nSeg);
128  }
129  }
130  }
131  }
132 
133  return 0;
134 }
135 
136 
137 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
Foam::polyLine
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:53
p
p
Definition: pEqn.H:62
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::IFstream
Input from file stream.
Definition: IFstream.H:81
Foam::argList::size
label size() const
Return the number of arguments.
Definition: argListI.H:84
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:97
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
CatmullRomSpline.H
Foam::BSpline::position
point position(const scalar lambda) const
Return the point position corresponding to the curve parameter.
Definition: BSpline.C:43
Foam::CatmullRomSpline
An implementation of Catmull-Rom splines (sometimes known as Overhauser splines).
Definition: CatmullRomSpline.H:75
Foam::BSpline
An implementation of B-splines.
Definition: BSpline.H:73
Foam::CatmullRomSpline::position
point position(const scalar lambda) const
Return the point position corresponding to the curve parameter.
Definition: CatmullRomSpline.C:43
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
Foam::argList::printUsage
void printUsage() const
Print usage.
Definition: argList.C:1045
IFstream.H
Foam::argList::optionLookupOrDefault
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:237
Foam::polyLine::position
point position(const scalar) const
Return the point position corresponding to the curve parameter.
Definition: polyLine.C:118
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
BSpline.H
printPoint
Ostream & printPoint(Ostream &os, const point &p)
Definition: Test-spline.C:36
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
vector.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
main
int main(int argc, char *argv[])
Definition: Test-spline.C:46
Foam::argList::noParallel
static void noParallel()
Remove the parallel options.
Definition: argList.C:161
args
Foam::argList args(argc, argv)
lambda
dimensionedScalar lambda(laminarTransport.lookup("lambda"))