Loading...
Searching...
No Matches
objective.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2020 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28Class
29 Foam::objective
30
31Description
32 Abstract base class for objective functions. No point in making this
33 runTime selectable since its children will have different constructors.
34
35SourceFiles
36 objective.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef objective_H
41#define objective_H
42
43#include "localIOdictionary.H"
44#include "autoPtr.H"
46#include "OFstream.H"
47#include "boundaryFieldsFwd.H"
48#include "solverControl.H"
49#include "objectiveFwd.H"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
56/*---------------------------------------------------------------------------*\
57 Class objective Declaration
58\*---------------------------------------------------------------------------*/
59
60class objective
61:
64protected:
66 // Protected data
68 const fvMesh& mesh_;
72 const word objectiveName_;
74 bool nullified_;
75 bool normalize_;
77
78 //- Objective function value and weight
79 scalar J_;
80
81 //- Average objective value
82 scalar JMean_;
83
84 //- Objective weight
85 scalar weight_;
87 //- Whether the objective is computed or not
88 // Some objective (e.g. geometric ones) might not change from one
89 // iteration of the primal solver to the next and might be expensive
90 // to evaluate. This can be used to compute them only once per
91 // optimisation cycle
92 bool computed_;
93
94 //- Normalization factor
97 //- Target value, in case the objective is used as a constraint
98 // Should be used with caution and with updateMethods
99 // than get affected by the target value, without
100 // requiring a sqr (e.g. SQP, MMA)
102
103 //- Target on the left hand-side of a double inequality,
104 //- for double sided constraints
106
107 //- Objective integration start and end times (for unsteady flows)
111 //- List of adjoint fields for which this objective will contribute
112 //- sources to their equations.
113 // Only for volume-based objectives for the moment
115
116 //- Contribution to field sensitivity derivatives
117 // Topology optimisation or other variants with
118 // as many design variables as the mesh cells
119 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
122 //- Contribution to sensitivity derivatives with a
123 //- random number of design variables
124 //- (neither surface, nor volume based)
125 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
126 autoPtr<scalarField> dJdbFieldPtr_;
127
128 // Contribution to surface sensitivity derivatives
129 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
131 //- Term from material derivative
133
134 //- Term multiplying delta(n dS)/delta b
136
137 //- Term multiplying delta(n)/delta b
139
140 //- Term multiplying delta(x)/delta b at the boundary
142
143 //- Term multiplying delta(x)/delta b at the boundary
144 //- for objectives that directly depend on x, e.g. moment
145 //- Needed in both FI and SI computations
147
148 //- Contribution located in specific parts of a patch.
149 //- Mainly intended for patch boundary edges contributions, e.g.
150 //- NURBS surfaces G1 continuity
152
153 // Contribution to volume-based sensitivities from volume-based
154 // objective functions
155 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
157 //- Multiplier of d(Volume)/db
159
160 //- Emerging from volume objectives that include spatial derivatives
162
163 //- Output file variables
165
166 //- File to keep the objective values after the end of the primal solver
168
169 //- File to keep the objective values at each iteration of the primal
170 //- solver
172
173 //- File to keep the average objective values after the end of the
174 //- primal solver
176
177 //- Default width of entries when writing in the objective files
178 unsigned int width_;
179
180
181 // Protected Member Functions
182
183 //- Set the output file ptr
184 void setObjectiveFilePtr() const;
186 //- Set the output file ptr for the instantaneous value
187 void setInstantValueFilePtr() const;
188
189 //- Set the output file ptr for the mean value
190 void setMeanValueFilePtr() const;
191
192
193private:
195 // Private Member Functions
196
197 //- No copy construct
198 objective(const objective&) = delete;
200 //- No copy assignment
201 void operator=(const objective&) = delete;
202
203 //- Make objective Function Folder
204 void makeFolder();
205
206
207public:
208
209 //- Runtime type information
210 TypeName("objective");
211
212
213 // Declare run-time constructor selection table
214
216 (
217 autoPtr,
218 objective,
219 objective,
220 (
221 const fvMesh& mesh,
222 const dictionary& dict,
223 const word& adjointSolverName,
224 const word& primalSolverName
225 ),
226 (mesh, dict, adjointSolverName, primalSolverName)
227 );
228
229
230 // Constructors
231
232 //- Construct from components
233 objective
234 (
235 const fvMesh& mesh,
236 const dictionary& dict,
237 const word& adjointSolverName,
238 const word& primalSolverName
239 );
240
241
242 // Selectors
243
244 //- Return a reference to the selected turbulence model
246 (
247 const fvMesh& mesh,
248 const dictionary& dict,
249 const word& objectiveType,
250 const word& adjointSolverName,
251 const word& primalSolverName
252 );
253
254
255 //- Destructor
256 virtual ~objective() = default;
257
258
259 // Member Functions
260
261 virtual bool readDict(const dictionary& dict);
262
263 //- Return the instantaneous objective function value
264 virtual scalar J() = 0;
265
266 //- Return the objective function of the optimisation cycle.
267 // This corresponds to the mean value, if it exists, or the
268 // instantaneous value otherwise
269 virtual scalar JCycle(bool negate = false) const;
270
271 //- Accumulate contribution for the mean objective value
272 // For steady-state runs
274
275
276 //- Accumulate contribution for the mean objective value
277 // For unsteady runs
278 virtual void accumulateJMean();
279
280 //- Return the objective function weight
281 virtual scalar weight() const;
282
283 //- Is the objective normalized
284 virtual bool normalize() const;
285
286 //- Normalize all fields allocated by the objective
287 virtual void doNormalization();
288
289 //- Check whether this is an objective integration time
290 virtual bool isWithinIntegrationTime() const;
291
292 //- Increment integration times
293 virtual void incrementIntegrationTimes(const scalar timeSpan);
294
295 //- Contribution to field sensitivities
296 inline const volScalarField& dJdb() const;
297
298 //- Contribution to sensitivities with a random number of designVars
299 inline const scalarField& dJdbField() const;
300
301 //- Contribution to surface sensitivities for a specific patch
302 inline const fvPatchVectorField& boundarydJdb(const label) const;
303
304 //- Multiplier of delta(n dS)/delta b
305 inline const fvPatchVectorField& dSdbMultiplier(const label) const;
306
307 //- Multiplier of delta(n dS)/delta b
308 inline const fvPatchVectorField& dndbMultiplier(const label) const;
309
310 //- Multiplier of delta(x)/delta b
311 inline const fvPatchVectorField& dxdbMultiplier(const label) const;
312
313 //- Multiplier of delta(x)/delta b
315 (
316 const label
317 ) const;
318
319 //- Multiplier located at patch boundary edges
321 (
322 const label patchI,
323 const label edgeI
324 ) const;
325
326 //- Contribution to surface sensitivities for all patches
327 inline const boundaryVectorField& boundarydJdb() const;
328
329 //- Multiplier of delta(n dS)/delta b for all patches
330 inline const boundaryVectorField& dSdbMultiplier() const;
331
332 //- Multiplier of delta(n dS)/delta b for all patches
333 inline const boundaryVectorField& dndbMultiplier() const;
335 //- Multiplier of delta(x)/delta b for all patches
336 inline const boundaryVectorField& dxdbMultiplier() const;
337
338 //- Multiplier of delta(x)/delta b for all patches
339 inline const boundaryVectorField& dxdbDirectMultiplier() const;
340
341 //- Multiplier located at patch boundary edges
342 inline const vectorField3& boundaryEdgeMultiplier() const;
343
344 //- Multiplier of grad( delta(x)/delta b) for volume-based sensitivities
345 inline const volScalarField& divDxDbMultiplier() const;
346
347 //- Multiplier of grad( delta(x)/delta b) for volume-based sensitivities
348 inline const volTensorField& gradDxDbMultiplier() const;
349
350 //- Update objective function derivatives
351 virtual void update() = 0;
352
353 //- Nullify adjoint contributions
354 virtual void nullify();
355
356 //- Update normalization factors, for objectives in
357 //- which the factor is not known a priori
358 virtual void updateNormalizationFactor();
359
360
361 virtual void update_dJdb()
362 {}
363
364 virtual void update_dJdbField()
365 {}
366
367 //- Update objective function derivative term
368 virtual void update_boundarydJdb()
369 {}
370
371 //- Update d (normal dS) / db multiplier. Surface-based sensitivity term
372 virtual void update_dSdbMultiplier()
373 {}
374
375 //- Update d (normal) / db multiplier. Surface-based sensitivity term
376 virtual void update_dndbMultiplier()
377 {}
378
379 //- Update d (x) / db multiplier. Surface-based sensitivity term
380 virtual void update_dxdbMultiplier()
381 {}
382
383 //- Update d (x) / db multiplier. Surface and volume-based sensitivity
384 //- term
385 virtual void update_dxdbDirectMultiplier()
386 {}
387
388 //- Update boundary edge contributions
390 {}
391
392 //- Update div( dx/db multiplier). Volume-based sensitivity term
393 virtual void update_divDxDbMultiplier()
394 {}
395
396 //- Update grad( dx/db multiplier). Volume-based sensitivity term
397 virtual void update_gradDxDbMultiplier()
398 {}
399
400
401 //- Manipulate fvVectorMatrix through the objective
402 virtual void addSource(fvVectorMatrix& matrix)
403 {}
404
405 //- Manipulate fvVectorMatrix through the objective
406 virtual void addSource(fvScalarMatrix& matrix)
407 {}
408
409 //- Write objective function history
410 virtual bool write(const bool valid = true) const;
411
412 //- Write objective function history at each primal solver iteration
413 virtual void writeInstantaneousValue() const;
414
415 //- Append a blank line after the end of optimisation cycle to the
416 //- file holding the instantaneous objective values to facilitate
417 //- plotting
418 virtual void writeInstantaneousSeparator() const;
419
420 //- Write mean objective function history
421 virtual void writeMeanValue() const;
422
423 //- Write averaged objective for continuation
424 virtual bool writeData(Ostream& os) const;
425
426 // Helper write functions
427
428 //- Write any information that needs to go the header of the file
429 // (e.g. targets, directions, etc)
430 virtual void addHeaderInfo() const;
431
432 //- Write headers for additional columns
433 virtual void addHeaderColumns() const;
434
435 //- Write information to additional columns
436 virtual void addColumnValues() const;
437
438 //- Return the objective name
439 inline const word& objectiveName() const;
440
441 //- Should the objective be written to file upon calling write()?
442 inline bool shouldWrite() const;
443
444 //- Should the objective be written to file upon calling write()?
445 inline void setWrite(const bool shouldWrite);
446
447 // Inline functions for checking whether pointers are set or not
448 inline bool hasdJdb() const noexcept;
449 inline bool hasdJdbField() const noexcept;
450 inline bool hasBoundarydJdb() const noexcept;
451 inline bool hasdSdbMult() const noexcept;
452 inline bool hasdndbMult() const noexcept;
453 inline bool hasdxdbMult() const noexcept;
454 inline bool hasdxdbDirectMult() const noexcept;
455 inline bool hasBoundaryEdgeContribution() const noexcept;
456 inline bool hasDivDxDbMult() const noexcept;
457 inline bool hasGradDxDbMult() const noexcept;
458
459 // Inline functions for checking whether integration times are set
460 inline bool hasIntegrationStartTime() const noexcept;
461 inline bool hasIntegrationEndTime() const noexcept;
462
463 //- Set the computed status of the objective
464 inline void setComputed(const bool isComputed) noexcept;
465
466 //- Return objective dictionary
467 const dictionary& dict() const;
468};
469
470
471// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
472
473} // End namespace Foam
474
475// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476
477#include "objectiveI.H"
478
479// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
480
481#endif
482
483// ************************************************************************* //
Useful typenames for fields defined only at the boundaries.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
virtual void writeInstantaneousSeparator() const
Append a blank line after the end of optimisation cycle to the file holding the instantaneous objecti...
Definition objective.C:602
bool hasdSdbMult() const noexcept
Definition objectiveI.H:181
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition objective.C:583
virtual void addSource(fvScalarMatrix &matrix)
Manipulate fvVectorMatrix through the objective.
Definition objective.H:553
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
Definition objective.H:515
autoPtr< OFstream > instantValueFilePtr_
File to keep the objective values at each iteration of the primal solver.
Definition objective.H:215
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
Definition objective.H:540
bool hasIntegrationStartTime() const noexcept
Definition objectiveI.H:223
autoPtr< boundaryVectorField > bdJdbPtr_
Contribution to field sensitivity derivatives.
Definition objective.H:156
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
Definition objective.H:522
void setComputed(const bool isComputed) noexcept
Set the computed status of the objective.
Definition objectiveI.H:235
virtual void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition objective.C:331
const fvMesh & mesh_
Definition objective.H:63
virtual void update()=0
Update objective function derivatives.
Definition objective.C:453
scalar JMean_
Average objective value.
Definition objective.H:81
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition objective.H:161
const boundaryVectorField & dxdbMultiplier() const
Multiplier of delta(x)/delta b for all patches.
Definition objectiveI.H:125
const word objectiveName_
Definition objective.H:67
virtual void addSource(fvVectorMatrix &matrix)
Manipulate fvVectorMatrix through the objective.
Definition objective.H:547
static autoPtr< objective > New(const fvMesh &mesh, const dictionary &dict, const word &objectiveType, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
Definition objective.C:221
virtual void update_dJdbField()
Definition objective.H:491
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition objective.C:77
virtual void writeMeanValue() const
Write mean objective function history.
Definition objective.C:614
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition objective.H:209
autoPtr< OFstream > meanValueFilePtr_
File to keep the average objective values after the end of the primal solver.
Definition objective.H:221
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
Definition objective.H:130
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Definition objective.H:534
virtual void addHeaderColumns() const
Write headers for additional columns.
Definition objective.C:659
bool hasGradDxDbMult() const noexcept
Definition objectiveI.H:217
bool hasBoundarydJdb() const noexcept
Definition objectiveI.H:175
bool hasdJdb() const noexcept
Definition objectiveI.H:163
virtual scalar weight() const
Return the objective function weight.
Definition objective.C:353
virtual bool writeData(Ostream &os) const
Write averaged objective for continuation.
Definition objective.C:641
autoPtr< scalar > targetLeft_
Target on the left hand-side of a double inequality, for double sided constraints.
Definition objective.H:116
const boundaryVectorField & dxdbDirectMultiplier() const
Multiplier of delta(x)/delta b for all patches.
Definition objectiveI.H:132
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x,...
Definition objective.H:178
virtual void update_boundarydJdb()
Update objective function derivative term.
Definition objective.H:497
declareRunTimeNewSelectionTable(autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
virtual void nullify()
Nullify adjoint contributions.
Definition objective.C:474
dictionary dict_
Definition objective.H:64
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows).
Definition objective.H:121
const word & objectiveName() const
Return the objective name.
Definition objectiveI.H:26
const scalarField & dJdbField() const
Contribution to sensitivities with a random number of designVars.
Definition objectiveI.H:50
bool hasIntegrationEndTime() const noexcept
Definition objectiveI.H:229
bool hasdndbMult() const noexcept
Definition objectiveI.H:187
void setObjectiveFilePtr() const
Set the output file ptr.
Definition objective.C:56
const boundaryVectorField & dSdbMultiplier() const
Multiplier of delta(n dS)/delta b for all patches.
Definition objectiveI.H:113
autoPtr< scalar > normFactor_
Normalization factor.
Definition objective.H:101
virtual scalar JCycle(bool negate=false) const
Return the objective function of the optimisation cycle.
Definition objective.C:264
scalar weight_
Objective weight.
Definition objective.H:86
bool computeMeanFields_
Definition objective.H:68
virtual bool readDict(const dictionary &dict)
Definition objective.C:257
unsigned int width_
Default width of entries when writing in the objective files.
Definition objective.H:226
virtual scalar J()=0
Return the instantaneous objective function value.
autoPtr< scalar > integrationEndTimePtr_
Definition objective.H:122
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
Definition objective.H:503
void setWrite(const bool shouldWrite)
Should the objective be written to file upon calling write()?
Definition objectiveI.H:38
const volScalarField & divDxDbMultiplier() const
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition objectiveI.H:151
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
Definition objective.C:653
bool hasdJdbField() const noexcept
Definition objectiveI.H:169
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
Definition objective.H:509
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition objective.H:194
const volScalarField & dJdb() const
Contribution to field sensitivities.
Definition objectiveI.H:44
bool hasDivDxDbMult() const noexcept
Definition objectiveI.H:211
scalar J_
Objective function value and weight.
Definition objective.H:76
virtual ~objective()=default
Destructor.
TypeName("objective")
Runtime type information.
bool hasBoundaryEdgeContribution() const noexcept
Definition objectiveI.H:205
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition objective.C:65
virtual void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition objective.C:437
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition objective.C:365
const boundaryVectorField & dndbMultiplier() const
Multiplier of delta(n dS)/delta b for all patches.
Definition objectiveI.H:119
virtual bool normalize() const
Is the objective normalized.
Definition objective.C:359
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
const word primalSolverName_
Definition objective.H:66
bool computed_
Whether the objective is computed or not.
Definition objective.H:96
virtual void updateNormalizationFactor()
Update normalization factors, for objectives in which the factor is not known a priori.
Definition objective.C:300
virtual void update_dJdb()
Definition objective.H:488
fileName objFunctionFolder_
Output file variables.
Definition objective.H:204
autoPtr< scalar > target_
Target value, in case the objective is used as a constraint.
Definition objective.H:110
const boundaryVectorField & boundarydJdb() const
Contribution to surface sensitivities for all patches.
Definition objectiveI.H:107
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition objective.H:166
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition objective.H:528
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition objective.H:171
const word adjointSolverName_
Definition objective.H:65
virtual void addColumnValues() const
Write information to additional columns.
Definition objective.C:665
bool hasdxdbDirectMult() const noexcept
Definition objectiveI.H:199
virtual bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition objective.C:415
bool shouldWrite() const
Should the objective be written to file upon calling write()?
Definition objectiveI.H:32
autoPtr< vectorField3 > bEdgeContribution_
Contribution located in specific parts of a patch. Mainly intended for patch boundary edges contribut...
Definition objective.H:185
bool hasdxdbMult() const noexcept
Definition objectiveI.H:193
const vectorField3 & boundaryEdgeMultiplier() const
Multiplier located at patch boundary edges.
Definition objectiveI.H:138
const volTensorField & gradDxDbMultiplier() const
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition objectiveI.H:157
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition objective.H:199
Base class for solver control classes.
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Field< Field< vectorField > > vectorField3
List< word > wordList
List of word.
Definition fileName.H:60
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
fvMatrix< vector > fvVectorMatrix
Field< vector > vectorField
Specialisation of Field<T> for vector.
const direction noexcept
Definition scalarImpl.H:265
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
volVectorField::Boundary boundaryVectorField
fvPatchField< vector > fvPatchVectorField
runTime write()
Macros to ease declaration of run-time selection tables.
#define declareRunTimeNewSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection for derived classes.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68