blockEdge.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "blockEdge.H"
27 #include "blockVertex.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 (
42  const pointField& points,
43  const label start,
44  const label end
45 )
46 :
47  points_(points),
48  start_(start),
49  end_(end)
50 {}
51 
52 
54 (
55  const dictionary& dict,
56  const label index,
57  const pointField& points,
58  Istream& is
59 )
60 :
61  points_(points),
62  start_(blockVertex::read(is, dict)),
63  end_(blockVertex::read(is, dict))
64 {}
65 
66 
68 {
70  return autoPtr<blockEdge>(nullptr);
71 }
72 
73 
75 (
76  const dictionary& dict,
77  const label index,
78  const searchableSurfaces& geometry,
79  const pointField& points,
80  Istream& is
81 )
82 {
83  if (debug)
84  {
85  InfoInFunction << "Constructing blockEdge" << endl;
86  }
87 
88  const word edgeType(is);
89 
90  IstreamConstructorTable::iterator cstrIter =
91  IstreamConstructorTablePtr_->find(edgeType);
92 
93  if (cstrIter == IstreamConstructorTablePtr_->end())
94  {
96  << "Unknown blockEdge type "
97  << edgeType << nl << nl
98  << "Valid blockEdge types are" << endl
99  << IstreamConstructorTablePtr_->sortedToc()
100  << abort(FatalError);
101  }
102 
103  return autoPtr<blockEdge>(cstrIter()(dict, index, geometry, points, is));
104 }
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 (
111  const pointField& points,
112  const label start,
113  const label end,
114  const pointField& otherKnots
115 )
116 {
117  pointField allKnots(otherKnots.size() + 2);
118 
119  // Start/end knots
120  allKnots[0] = points[start];
121  allKnots[otherKnots.size() + 1] = points[end];
122 
123  // Intermediate knots
124  forAll(otherKnots, knotI)
125  {
126  allKnots[knotI+1] = otherKnots[knotI];
127  }
128 
129  return allKnots;
130 }
131 
132 
135 {
136  tmp<pointField> tpoints(new pointField(lambdas.size()));
137  pointField& points = tpoints.ref();
138 
139  forAll(lambdas, i)
140  {
141  points[i] = position(lambdas[i]);
142  }
143  return tpoints;
144 }
145 
146 
147 void Foam::blockEdge::write(Ostream& os, const dictionary& d) const
148 {
149  blockVertex::write(os, start_, d);
150  os << tab;
151  blockVertex::write(os, end_, d);
152  os << endl;
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157 
159 {
160  os << p.start_ << tab << p.end_ << endl;
161 
162  return os;
163 }
164 
165 
166 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Define a curved edge that is parameterised for 0<lambda<1 between the start and end point.
Definition: blockEdge.H:57
static autoPtr< blockEdge > New(const dictionary &dict, const label index, const searchableSurfaces &geometry, const pointField &, Istream &)
New function which constructs and returns pointer to a blockEdge.
Definition: blockEdge.C:75
blockEdge(const pointField &points, const label start, const label end)
Construct from components.
Definition: blockEdge.C:41
virtual autoPtr< blockEdge > clone() const
Clone function.
Definition: blockEdge.C:67
virtual point position(const scalar) const =0
Return the point position corresponding to the curve parameter.
void write(Ostream &, const dictionary &) const
Write edge with variable backsubstitution.
Definition: blockEdge.C:147
static pointField appendEndPoints(const pointField &, const label start, const label end, const pointField &otherKnots)
Return a complete point field by appending the start/end points.
Definition: blockEdge.C:110
Define a block vertex.
Definition: blockVertex.H:49
static void write(Ostream &, const label, const dictionary &)
Write vertex index with optional name backsubstitution.
Definition: blockVertex.C:120
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Container for searchableSurfaces.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
static const char tab
Definition: Ostream.H:265
errorManip< error > abort(error &err)
Definition: errorManip.H:131
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
defineTypeNameAndDebug(combustionModel, 0)
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p