Loading...
Searching...
No Matches
sensitivitySurfacePoints.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2020, 2022 PCOpt/NTUA
9 Copyright (C) 2013-2020, 2022 FOSS GP
10 Copyright (C) 2019-2022 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
28\*---------------------------------------------------------------------------*/
29
31#include "deltaBoundary.H"
32#include "designVariables.H"
33#include "syncTools.H"
34#include "symmetryFvPatch.H"
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40namespace Foam
41{
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
47(
49 sensitivitySurfacePoints,
51);
52
53// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54
56{
57 // Populate extendedPatchIDs
58 label pI(0);
59 labelList extendedPatchIDs(mesh_.boundary().size(), -1);
60 forAll(mesh_.boundary(), patchI)
61 {
62 const fvPatch& pp = mesh_.boundary()[patchI];
63 bool isSymmetry
65 if (!isA<coupledFvPatch>(pp) && !isA<emptyFvPatch>(pp) && !isSymmetry)
66 {
67 extendedPatchIDs[pI++] = patchI;
68 }
69 }
70 extendedPatchIDs.setSize(pI);
71 return labelHashSet(extendedPatchIDs);
72}
73
74
77 word suffix(adjointMeshMovementSolver_ ? "ESI" : "SI");
78 suffix = suffix + word(dict().getOrDefault<word>("suffix", word::null));
79 setSuffix(adjointSolver_.solverName() + suffix);
80}
81
82
84{
85 // List with mesh faces. Global addressing
86 const faceList& faces = mesh_.faces();
87
88 // Geometry differentiation engine
89 deltaBoundary dBoundary(mesh_);
90
91 for (const label patchI : extendedPatchIDs_)
92 {
93 const fvPatch& patch = mesh_.boundary()[patchI];
94 vectorField nf(patch.nf());
95
96 // Point sens result for patch
97 vectorField& pointPatchSens = wallPointSensVecPtr_()[patchI];
98
99 // Face sens for patch
100 vectorField facePatchSens = dxdbMult_()[patchI];
101 if (dxdbDirectMult_)
102 {
103 facePatchSens += dxdbDirectMult_()[patchI];
104 }
105 if (bcDxDbMult_)
106 {
107 facePatchSens += bcDxDbMult_()[patchI];
108 }
109
110 // Correspondance of local point addressing to global point addressing
111 const labelList& meshPoints = patch.patch().meshPoints();
112
113 // Each local patch point belongs to these local patch faces
114 // (local numbering)
115 const labelListList& patchPointFaces = patch.patch().pointFaces();
116
117 // Index of first face in patch
118 const label patchStartIndex = patch.start();
119
120 // Loop over patch points.
121 // Collect contributions from each boundary face this point belongs to
122 forAll(meshPoints, ppI)
123 {
124 const labelList& pointFaces = patchPointFaces[ppI];
125 forAll(pointFaces, pfI)
126 {
127 label localFaceIndex = pointFaces[pfI];
128 label globalFaceIndex = patchStartIndex + localFaceIndex;
129 const face& faceI = faces[globalFaceIndex];
130
131 // Point coordinates. All indices in global numbering
132 pointField p(faceI.points(mesh_.points()));
133 tensorField p_d(faceI.size(), Zero);
134 forAll(faceI, facePointI)
135 {
136 if (faceI[facePointI] == meshPoints[ppI])
137 {
138 p_d[facePointI] = tensor::I;
139 }
140 }
141 tensorField deltaNormals =
142 dBoundary.makeFaceCentresAndAreas_d(p, p_d);
143
144 if (isSymmetryPoint_[meshPoints[ppI]])
145 {
146 const vector& n = symmPointNormal_[meshPoints[ppI]];
147 deltaNormals =
148 //0.5*(deltaNormals + transform(I - 2.0*sqr(n), deltaNormals));
149 (deltaNormals + transform(I - 2.0*sqr(n), deltaNormals));
150 }
151
152 // Element [0] is the variation in the face center
153 // (dxFace/dxPoint)
154 const tensor& deltaCf = deltaNormals[0];
155 pointPatchSens[ppI] += facePatchSens[localFaceIndex] & deltaCf;
156
157 // Term multiplying d(Sf)/d(point displacement) and
158 // d(nf)/d(point displacement)
159 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
160 // Element [1] is the variation in the (dimensional) normal
161 if (dSfdbMult_)
162 {
163 const tensor& deltaSf = deltaNormals[1];
164 pointPatchSens[ppI] +=
165 dSfdbMult_()[patchI][localFaceIndex] & deltaSf;
166 }
167
168 // Element [2] is the variation in the unit normal
169 if (dnfdbMult_)
170 {
171 const tensor& deltaNf = deltaNormals[2];
172 pointPatchSens[ppI] +=
173 dnfdbMult_()[patchI][localFaceIndex] & deltaNf;
175 }
176 }
177 }
178}
179
180
182(
183 vectorField& pointNormals,
184 scalarField& pointMagSf
185)
186{
187 for (const label patchI : extendedPatchIDs_)
188 {
189 const fvPatch& patch = mesh_.boundary()[patchI];
190 const scalarField& magSf = patch.magSf();
191 vectorField nf(patch.nf());
192
193 // Correspondance of local point addressing to global point addressing
194 const labelList& meshPoints = patch.patch().meshPoints();
195
196 // Each local patch point belongs to these local patch faces
197 // (local numbering)
198 const labelListList& patchPointFaces = patch.patch().pointFaces();
199
200 // Loop over patch points
201 forAll(meshPoints, ppI)
202 {
203 const labelList& pointFaces = patchPointFaces[ppI];
204 forAll(pointFaces, pfI)
205 {
206 const label localFaceIndex = pointFaces[pfI];
207
208 // Accumulate information for point normals
209 pointNormals[meshPoints[ppI]] += nf[localFaceIndex];
210 pointMagSf[meshPoints[ppI]] += magSf[localFaceIndex];
211 }
212 }
213 }
214
216 (
217 mesh_,
218 pointNormals,
219 plusEqOp<vector>(),
221 );
223 (
224 mesh_,
225 pointMagSf,
226 plusEqOp<scalar>(),
227 scalar(0)
228 );
229
231 {
232 pointScalarField MagSf
233 (
234 IOobject
235 (
236 "pointMagSf",
237 mesh_.time().timeName(),
238 mesh_,
241 ),
244 );
246 (
247 IOobject
248 (
249 "pointNf",
250 mesh_.time().timeName(),
251 mesh_,
254 ),
257 );
258 MagSf.primitiveFieldRef() = pointMagSf;
259 Nf.primitiveFieldRef() = pointNormals;
260 Nf.primitiveFieldRef().normalise();
261 MagSf.write();
262 Nf.write();
263 }
264}
265
266
268{
269 // Allocate appropriate space for sensitivities
270 label nTotalPoints(0);
271 for (const label patchI : sensitivityPatchIDs_)
272 {
273 nTotalPoints += mesh_.boundaryMesh()[patchI].nPoints();
274 }
275
276 // Derivatives for all (x,y,z) components of the displacement
277 derivatives_ = scalarField(3*nTotalPoints, Zero);
278}
279
280
281// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
282
283sensitivitySurfacePoints::sensitivitySurfacePoints
284(
285 const fvMesh& mesh,
286 const dictionary& dict,
288)
289:
291 writeGeometricInfo_(false),
292 includeSurfaceArea_(false),
293 isSymmetryPoint_(mesh.nPoints(), false),
294 symmPointNormal_(mesh.nPoints(), Zero),
295 extendedPatchIDs_(populateExtendedIDs())
296{
297 if (debug)
298 {
299 Info<< "Extended sensitivity patches " << nl;
300 for (const label patchI : extendedPatchIDs_)
301 {
302 Info<< mesh_.boundary()[patchI].name() << endl;
303 }
304 }
305 read();
307
308 // Allocate boundary field pointer
311 (
313 );
315 (
317 );
318
320
321 // Populate symmetry patches
322 forAll(mesh_.boundary(), patchI)
323 {
324 const fvPatch& pp = mesh_.boundary()[patchI];
325 bool isSymmetry
327 if (isSymmetry)
328 {
329 const labelList& meshPoints = pp.patch().meshPoints();
330 const vectorField& pointNormals = pp.patch().pointNormals();
331 forAll(meshPoints, pI)
332 {
333 const label pointi = meshPoints[pI];
334 isSymmetryPoint_[pointi] = true;
335 symmPointNormal_[pointi] = pointNormals[pI];
336 }
338 }
339}
340
341
342// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
343
345{
347 dict().getOrDefault<bool>("writeGeometricInfo", false);
348 // Point sensitivities do not include the surface area by default
350 dict().getOrDefault<bool>("includeSurfaceArea", false);
351}
352
353
355{
357 {
358 read();
359 return true;
361
362 return false;
363}
364
365
368{
369 return extendedPatchIDs_;
370}
371
372
374(
375 autoPtr<designVariables>& designVars
376)
377{
378 // Make sure we have the proper size for the sensitivities
380
381 // Assemble the multipliers of dxdbFace, as in the ESI approach
383
384 // Geometric (or "direct") sensitivities are better computed directly on
385 // the points. Compute them and add the ones that depend on dxFace/dxPoint
387
388 // polyPatch::pointNormals will give the wrong result for points
389 // belonging to multiple patches or patch-processorPatch intersections.
390 // Keeping a mesh-wide field to allow easy reduction using syncTools.
391 // A bit expensive? Better way?
392 vectorField pointNormals(mesh_.nPoints(), Zero);
393 scalarField pointMagSf(mesh_.nPoints(), Zero);
394 constructGlobalPointNormalsAndAreas(pointNormals, pointMagSf);
395
396 // Do parallel communications to avoid wrong values at processor boundaries
397 // Global field for accumulation
398 vectorField pointSensGlobal(mesh_.nPoints(), Zero);
399 for (const label patchI : extendedPatchIDs_)
400 {
401 const labelList& meshPoints = mesh_.boundaryMesh()[patchI].meshPoints();
402 forAll(meshPoints, ppI)
403 {
404 const label globaPointI = meshPoints[ppI];
405 pointSensGlobal[globaPointI] += wallPointSensVecPtr_()[patchI][ppI];
406 }
407 }
408
409 /*
410 // Remove components normal to symmetry planes
411 forAll(mesh_.boundary(), patchI)
412 {
413 const fvPatch& patch = mesh_.boundary()[patchI];
414 if (isA<symmetryFvPatch>(patch) || isA<symmetryPlaneFvPatch>(patch))
415 {
416 // Deliberately using local point normals instead of the global ones,
417 // to get the direction normal to the symmetry plane itself
418 const vectorField& pn = patch.patch().pointNormals();
419 const labelList& meshPoints = patch.patch().meshPoints();
420 forAll(meshPoints, pI)
421 {
422 const label gpI = meshPoints[pI];
423 pointSensGlobal[gpI] -= wallPointSensVecPtr_()[patchI][pI];
424 pointSensGlobal[gpI] -=
425 (pointSensGlobal[gpI] & pn[pI])*pn[pI];
426 pointSensGlobal[gpI] *= 2;
427 }
428 }
429 }
430 */
431
432 // Accumulate dJ/dx_i
434 (
435 mesh_,
436 pointSensGlobal,
437 plusEqOp<vector>(),
439 );
440
441 // Transfer back to local fields
442 for (const label patchI : extendedPatchIDs_)
443 {
444 const labelList& meshPoints =
445 mesh_.boundaryMesh()[patchI].meshPoints();
446 wallPointSensVecPtr_()[patchI].map(pointSensGlobal, meshPoints);
447 }
448
449 // Compute normal sens and append to return field
450 label nPassedDVs(0);
451 const Vector<label>& sd = mesh_.solutionD();
452 for (const label patchI : sensitivityPatchIDs_)
453 {
454 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
455 //if (patch.size()>0)
456 {
457 const labelList& meshPoints = patch.meshPoints();
458
459 // Avoid storing unit point normals in the global list since we
460 // might divide multiple times with the number of faces belonging
461 // to the point. Instead do the division locally, per patch use
462 vectorField patchPointNormals(pointNormals, meshPoints);
463 patchPointNormals.normalise();
465 {
466 wallPointSensVecPtr_()[patchI] /=
467 scalarField(pointMagSf, meshPoints);
468 }
469 wallPointSensNormalPtr_()[patchI] =
470 wallPointSensVecPtr_()[patchI] & patchPointNormals;
472 wallPointSensNormalPtr_()[patchI] *patchPointNormals;
473
474 forAll(patch.localPoints(), pi)
475 {
476 const label gpi = nPassedDVs + pi;
477 const vector& pSens = wallPointSensVecPtr_()[patchI][pi];
478 derivatives_[3*gpi ] = scalar(sd[0] == -1 ? 0 : pSens.x());
479 derivatives_[3*gpi + 1] = scalar(sd[1] == -1 ? 0 : pSens.y());
480 derivatives_[3*gpi + 2] = scalar(sd[2] == -1 ? 0 : pSens.z());
481 }
482 nPassedDVs += patch.nPoints();
483 }
484 }
485
486 // Write derivative fields
487 write();
488
489 // Get processed sensitivities from designVariables, if present
490 if (designVars)
491 {
493 }
494}
495
496
497void sensitivitySurfacePoints::write(const word& baseName)
498{
501
503 {
504 volVectorField nfOnPatch
505 (
506 IOobject
507 (
508 "nfOnPatch",
509 mesh_.time().timeName(),
510 mesh_,
513 ),
514 mesh_,
515 Zero
516 );
517
518 volVectorField SfOnPatch
519 (
520 IOobject
521 (
522 "SfOnPatch",
523 mesh_.time().timeName(),
524 mesh_,
527 ),
528 mesh_,
529 Zero
530 );
531
532 volVectorField CfOnPatch
533 (
534 IOobject
535 (
536 "CfOnPatch",
537 mesh_.time().timeName(),
538 mesh_,
541 ),
542 mesh_,
543 Zero
544 );
545 for (const label patchI : sensitivityPatchIDs_)
546 {
547 const fvPatch& patch = mesh_.boundary()[patchI];
548 nfOnPatch.boundaryFieldRef()[patchI] = patch.nf();
549 SfOnPatch.boundaryFieldRef()[patchI] = patch.Sf();
550 CfOnPatch.boundaryFieldRef()[patchI] = patch.Cf();
551 }
552 }
553}
554
555
556// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
557
558} // End namespace Foam
559
560// ************************************************************************* //
constexpr scalar pi(M_PI)
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
Definition Field.C:600
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
autoPtr< pointBoundaryVectorField > wallPointSensNormalVecPtr_
Normal sens as vectors.
labelHashSet sensitivityPatchIDs_
Patches on which to compute shape sensitivities.
autoPtr< pointBoundaryVectorField > wallPointSensVecPtr_
Wall point sens w.r.t. (x,y.z).
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
autoPtr< pointBoundaryScalarField > wallPointSensNormalPtr_
Wall point sens projected to normal.
static const Tensor I
Definition Tensor.H:81
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
Abstract base class for adjoint-based sensitivities.
autoPtr< boundaryVectorField > dxdbDirectMult_
Multiplier of dCf/db, found in the objective function.
void setSuffix(const word &suffix)
Set suffix.
autoPtr< boundaryVectorField > bcDxDbMult_
Multiplier of dx/db, coming from boundary conditions that depend on the geometry, like rotatingWallVe...
adjointSolver & adjointSolver_
Reference to the underlaying adjoint solver.
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
autoPtr< boundaryVectorField > dxdbMult_
Multiplier of face dx/db.
scalarField derivatives_
The sensitivity derivative values.
autoPtr< boundaryVectorField > dnfdbMult_
Multiplier of dnf/db.
autoPtr< adjointMeshMovementSolver > adjointMeshMovementSolver_
Adjoint grid displacement solver.
autoPtr< boundaryVectorField > dSfdbMult_
Multiplier of dSf/db.
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Differentiation of the mesh data structure.
vectorField makeFaceCentresAndAreas_d(const pointField &p, const pointField &p_d)
Given a face and the points to be moved in the normal direction, find faceArea, faceCentre and unitVe...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition faceI.H:80
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nPoints() const noexcept
Number of mesh points.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Class for computing sensitivity derivatives using the Enhanced Surface Integral (E-SI) formulation,...
void computeDxDbMult()
Compute dxdbMult from its various components.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
Calculation of adjoint-based sensitivities at wall points using the E-SI formulation.
void computePointDerivativesSize()
Allocate the proper size for the point-based sensitivities.
void constructGlobalPointNormalsAndAreas(vectorField &pointNormals, scalarField &pointMagSf)
Construct globally correct point normals and point areas.
void setSuffixName()
Set suffix name for sensitivity fields.
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
bool includeSurfaceArea_
Include surface area in sens computation.
bool writeGeometricInfo_
Write geometric info for use by external programs.
void read()
Read controls and update solver pointers if necessary.
boolList isSymmetryPoint_
Is point belonging to a symmetry{Plane}.
labelHashSet extendedPatchIDs_
Extended patchIDs.
labelHashSet populateExtendedIDs() const
Set suffix name for sensitivity fields.
vectorField symmPointNormal_
Local point normal per symmetry point.
void finalisePointSensitivities()
Converts face sensitivities to point sensitivities and adds the ones directly computed in points (i....
virtual const labelHashSet & geometryVariationIntegrationPatches() const
Return set of patches on which to compute direct sensitivities.
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
const fvMesh & mesh_
Definition sensitivity.H:68
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
label nPoints
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
static const Identity< scalar > I
Definition Identity.H:100
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
autoPtr< List< Field< Type > > > createZeroBoundaryPointFieldPtr(const fvMesh &mesh, bool printAllocation=false)
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
runTime write()
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299