sampledSurfaces.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "sampledSurfaces.H"
30#include "polySurface.H"
31
32#include "mapPolyMesh.H"
33#include "volFields.H"
34#include "surfaceFields.H"
35#include "HashOps.H"
36#include "ListOps.H"
37#include "Time.H"
38#include "IndirectList.H"
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46
48 (
52 );
53}
54
55Foam::scalar Foam::sampledSurfaces::mergeTol_ = 1e-10;
56
57
58// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59
60Foam::polySurface* Foam::sampledSurfaces::getRegistrySurface
61(
62 const sampledSurface& s
63) const
64{
65 return s.getRegistrySurface
66 (
68 IOobject::groupName(name(), s.name())
69 );
70}
71
72
73Foam::polySurface* Foam::sampledSurfaces::storeRegistrySurface
74(
75 const sampledSurface& s
76)
77{
78 return s.storeRegistrySurface
79 (
80 storedObjects(),
81 IOobject::groupName(name(), s.name())
82 );
83}
84
85
86bool Foam::sampledSurfaces::removeRegistrySurface
87(
88 const sampledSurface& s
89)
90{
91 return s.removeRegistrySurface
92 (
93 storedObjects(),
94 IOobject::groupName(name(), s.name())
95 );
96}
97
98
99Foam::IOobjectList Foam::sampledSurfaces::preCheckFields()
100{
101 wordList allFields; // Just needed for warnings
102 HashTable<wordHashSet> selected;
103
104 IOobjectList objects(0);
105
106 if (loadFromFiles_)
107 {
108 // Check files for a particular time
109 objects = IOobjectList(obr_, obr_.time().timeName());
110
111 allFields = objects.names();
112 selected = objects.classes(fieldSelection_);
113 }
114 else
115 {
116 // Check currently available fields
117 allFields = obr_.names();
118 selected = obr_.classes(fieldSelection_);
119 }
120
121 // Parallel consistency (no-op in serial)
122 Pstream::mapCombineAllGather(selected, HashSetOps::plusEqOp<word>());
123
124
125 DynamicList<label> missed(fieldSelection_.size());
126
127 // Detect missing fields
128 forAll(fieldSelection_, i)
129 {
130 if (!ListOps::found(allFields, fieldSelection_[i]))
131 {
132 missed.append(i);
133 }
134 }
135
136 if (missed.size())
137 {
139 << nl
140 << "Cannot find "
141 << (loadFromFiles_ ? "field file" : "registered field")
142 << " matching "
143 << UIndirectList<wordRe>(fieldSelection_, missed) << endl;
144 }
145
146
147 // Currently only support volume and surface field types
148 label nVolumeFields = 0;
149 label nSurfaceFields = 0;
150
151 forAllConstIters(selected, iter)
152 {
153 const word& clsName = iter.key();
154 const label n = iter.val().size();
155
156 if (fieldTypes::volume.found(clsName))
157 {
158 nVolumeFields += n;
159 }
160 else if (fieldTypes::surface.found(clsName))
161 {
162 nSurfaceFields += n;
163 }
164 }
165
166 // Now propagate field counts (per surface)
167 // - can update writer even when not writing without problem
168
169 forAll(*this, surfi)
170 {
171 const sampledSurface& s = (*this)[surfi];
172 surfaceWriter& outWriter = writers_[surfi];
173
174 outWriter.nFields
175 (
176 nVolumeFields
177 + (s.withSurfaceFields() ? nSurfaceFields : 0)
178 +
179 (
180 // Face-id field, but not if the writer manages that itself
181 !s.isPointData() && s.hasFaceIds() && !outWriter.usesFaceIds()
182 ? 1 : 0
183 )
184 );
185 }
186
187 return objects;
188}
189
190
191Foam::autoPtr<Foam::surfaceWriter> Foam::sampledSurfaces::newWriter
192(
193 word writeType,
194 const dictionary& formatOptions,
195 const dictionary& surfDict
196)
197{
198 // Per-surface adjustment
199 surfDict.readIfPresent<word>("surfaceFormat", writeType);
200
201 dictionary options = formatOptions.subOrEmptyDict(writeType);
202
203 options.merge
204 (
205 surfDict.subOrEmptyDict("formatOptions").subOrEmptyDict(writeType)
206 );
207
208 return surfaceWriter::New(writeType, options);
209}
210
211
212// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
213
215(
216 const word& name,
217 const Time& runTime,
218 const dictionary& dict
219)
220:
221 functionObjects::fvMeshFunctionObject(name, runTime, dict),
223 loadFromFiles_(false),
224 verbose_(false),
225 onExecute_(false),
226 outputPath_
227 (
228 time_.globalPath()/functionObject::outputPrefix/name
229 ),
230 fieldSelection_(),
231 sampleFaceScheme_(),
232 sampleNodeScheme_(),
233 writers_(),
234 actions_(),
235 nFaces_()
236{
237 outputPath_.clean(); // Remove unneeded ".."
238
239 read(dict);
240}
241
242
244(
245 const word& name,
246 const objectRegistry& obr,
247 const dictionary& dict,
248 const bool loadFromFiles
249)
250:
251 functionObjects::fvMeshFunctionObject(name, obr, dict),
253 loadFromFiles_(loadFromFiles),
254 verbose_(false),
255 onExecute_(false),
256 outputPath_
257 (
258 time_.globalPath()/functionObject::outputPrefix/name
259 ),
260 fieldSelection_(),
261 sampleFaceScheme_(),
262 sampleNodeScheme_(),
263 writers_(),
264 actions_(),
265 nFaces_()
266{
267 outputPath_.clean(); // Remove unneeded ".."
268
269 read(dict);
270}
271
272
273// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
274
275bool Foam::sampledSurfaces::verbose(const bool on) noexcept
276{
277 bool old(verbose_);
278 verbose_ = on;
279 return old;
280}
281
282
284{
285 fvMeshFunctionObject::read(dict);
286
288 writers_.clear();
289 actions_.clear();
290 nFaces_.clear();
291 fieldSelection_.clear();
292
293 verbose_ = dict.getOrDefault("verbose", false);
294 onExecute_ = dict.getOrDefault("sampleOnExecute", false);
295
296 sampleFaceScheme_ =
297 dict.getOrDefault<word>("sampleScheme", "cell");
298
299 sampleNodeScheme_ =
300 dict.getOrDefault<word>("interpolationScheme", "cellPoint");
301
302 const entry* eptr = dict.findEntry("surfaces");
303
304 // Surface writer type and format options
305 const word writerType =
306 (eptr ? dict.get<word>("surfaceFormat") : word::null);
307
308 const dictionary formatOptions(dict.subOrEmptyDict("formatOptions"));
309
310 // Store on registry?
311 const bool dfltStore = dict.getOrDefault("store", false);
312
313 if (eptr && eptr->isDict())
314 {
315 PtrList<sampledSurface> surfs(eptr->dict().size());
316
317 actions_.resize(surfs.size(), ACTION_WRITE); // Default action
318 writers_.resize(surfs.size());
319 nFaces_.resize(surfs.size(), Zero);
320
321 label surfi = 0;
322
323 for (const entry& dEntry : eptr->dict())
324 {
325 if (!dEntry.isDict())
326 {
327 continue;
328 }
329
330 const dictionary& surfDict = dEntry.dict();
331
334 (
335 dEntry.keyword(),
336 mesh_,
337 surfDict
338 );
339
340 if (!surf || !surf->enabled())
341 {
342 continue;
343 }
344
345 // Define the surface
346 surfs.set(surfi, surf);
347
348 // Define additional action(s)
349 if (surfDict.getOrDefault("store", dfltStore))
350 {
351 actions_[surfi] |= ACTION_STORE;
352 }
353
354 // Define surface writer, but do NOT yet attach a surface
355 writers_.set
356 (
357 surfi,
358 newWriter(writerType, formatOptions, surfDict)
359 );
360
361 writers_[surfi].isPointData(surfs[surfi].isPointData());
362
363 // Use outputDir/TIME/surface-name
364 writers_[surfi].useTimeDir(true);
365 writers_[surfi].verbose(verbose_);
366
367 ++surfi;
368 }
369
370 surfs.resize(surfi);
371 actions_.resize(surfi);
372 writers_.resize(surfi);
373 surfaces().transfer(surfs);
374 }
375 else if (eptr)
376 {
377 // This is slightly trickier.
378 // We want access to the individual dictionaries used for construction
379
381
383 (
384 eptr->stream(),
385 sampledSurface::iNewCapture(mesh_, capture)
386 );
387
388 PtrList<sampledSurface> surfs(input.size());
389
390 actions_.resize(surfs.size(), ACTION_WRITE); // Default action
391 writers_.resize(surfs.size());
392 nFaces_.resize(surfs.size(), Zero);
393
394 label surfi = 0;
395
396 forAll(input, inputi)
397 {
398 const dictionary& surfDict = capture[inputi];
399
400 autoPtr<sampledSurface> surf = input.release(inputi);
401
402 if (!surf || !surf->enabled())
403 {
404 continue;
405 }
406
407 // Define the surface
408 surfs.set(surfi, surf);
409
410 // Define additional action(s)
411 if (surfDict.getOrDefault("store", dfltStore))
412 {
413 actions_[surfi] |= ACTION_STORE;
414 }
415
416 // Define surface writer, but do NOT yet attach a surface
417 writers_.set
418 (
419 surfi,
420 newWriter(writerType, formatOptions, surfDict)
421 );
422
423 writers_[surfi].isPointData(surfs[surfi].isPointData());
424
425 // Use outputDir/TIME/surface-name
426 writers_[surfi].useTimeDir(true);
427 writers_[surfi].verbose(verbose_);
428
429 ++surfi;
430 }
431
432 surfs.resize(surfi);
433 actions_.resize(surfi);
434 writers_.resize(surfi);
435 surfaces().transfer(surfs);
436 }
437
438
439 const auto& surfs = surfaces();
440
441 // Have some surfaces, so sort out which fields are needed and report
442
443 if (surfs.size())
444 {
445 nFaces_.resize(surfs.size(), Zero);
446
447 dict.readEntry("fields", fieldSelection_);
448 fieldSelection_.uniq();
449
450 forAll(*this, surfi)
451 {
452 const sampledSurface& s = (*this)[surfi];
453
454 if (!surfi)
455 {
456 Info<< "Sampled surface:" << nl;
457 }
458
459 Info<< " " << s.name() << " -> " << writers_[surfi].type();
460 if (actions_[surfi] & ACTION_STORE)
461 {
462 Info<< ", store on registry ("
463 << IOobject::groupName(name(), s.name()) << ')';
464 }
465 Info<< nl;
466 Info<< " ";
467 s.print(Info, 0);
468 Info<< nl;
469 }
470 Info<< nl;
471 }
472
473 if (debug && Pstream::master())
474 {
475 Pout<< "sample fields:" << fieldSelection_ << nl
476 << "sample surfaces:" << nl << '(' << nl;
477
478 for (const sampledSurface& s : surfaces())
479 {
480 Pout<< " " << s << nl;
481 }
482 Pout<< ')' << endl;
483 }
484
485 // Ensure all surfaces and merge information are expired
486 expire(true);
487
488 return true;
489}
490
491
492bool Foam::sampledSurfaces::performAction(unsigned request)
493{
494 // Update surfaces and store
495 bool ok = false;
496
497 forAll(*this, surfi)
498 {
499 sampledSurface& s = (*this)[surfi];
500
501 if (request & actions_[surfi])
502 {
503 if (s.update())
504 {
505 writers_[surfi].expire();
506 }
507
508 nFaces_[surfi] = returnReduce(s.faces().size(), sumOp<label>());
509
510 ok = ok || nFaces_[surfi];
511
512
513 // Store surfaces (even empty ones) otherwise we miss geometry
514 // updates.
515 // Any associated fields will be removed if the size changes
516
517 if ((request & actions_[surfi]) & ACTION_STORE)
518 {
519 storeRegistrySurface(s);
520 }
521 }
522 }
523
524 if (!ok)
525 {
526 // No surface with an applicable action or with faces to sample
527 return true;
528 }
529
530
531 // Determine availability of fields.
532 // Count per-surface number of fields, including Ids etc
533 // which only seems to be needed for VTK legacy
534
535 IOobjectList objects = preCheckFields();
536
537 // Update writers
538
539 forAll(*this, surfi)
540 {
541 const sampledSurface& s = (*this)[surfi];
542
543 if (((request & actions_[surfi]) & ACTION_WRITE) && nFaces_[surfi])
544 {
545 surfaceWriter& outWriter = writers_[surfi];
546
547 if (outWriter.needsUpdate())
548 {
549 outWriter.setSurface(s);
550 }
551
552 outWriter.open(outputPath_/s.name());
553
554 outWriter.beginTime(obr_.time());
555
556 // Face-id field, but not if the writer manages that itself
557 if (!s.isPointData() && s.hasFaceIds() && !outWriter.usesFaceIds())
558 {
559 // This is still quite horrible.
560
561 Field<label> ids(s.faceIds());
562
563 if
564 (
566 (
567 !ListOps::found(ids, lessOp1<label>(0)),
568 andOp<bool>()
569 )
570 )
571 {
572 // From 0-based to 1-based, provided there are no
573 // negative ids (eg, encoded solid/side)
574
575 ids += label(1);
576 }
577
578 writeSurface(outWriter, ids, "Ids");
579 }
580 }
581 }
582
583 // Sample fields
584
585 performAction<volScalarField>(objects, request);
586 performAction<volVectorField>(objects, request);
587 performAction<volSphericalTensorField>(objects, request);
588 performAction<volSymmTensorField>(objects, request);
589 performAction<volTensorField>(objects, request);
590
591 // Only bother with surface fields if a sampler supports them
592 if
593 (
594 testAny
595 (
596 surfaces(),
597 [] (const sampledSurface& s) { return s.withSurfaceFields(); }
598 )
599 )
600 {
601 performAction<surfaceScalarField>(objects, request);
602 performAction<surfaceVectorField>(objects, request);
603 performAction<surfaceSphericalTensorField>(objects, request);
604 performAction<surfaceSymmTensorField>(objects, request);
605 performAction<surfaceTensorField>(objects, request);
606 }
607
608
609 // Finish this time step
610 forAll(writers_, surfi)
611 {
612 if (((request & actions_[surfi]) & ACTION_WRITE) && nFaces_[surfi])
613 {
614 // Write geometry if no fields were written so that we still
615 // can have something to look at
616
617 if (!writers_[surfi].wroteData())
618 {
619 writers_[surfi].write();
620 }
621
622 writers_[surfi].endTime();
623 }
624 }
625
626 return true;
627}
628
629
631{
632 if (onExecute_)
633 {
634 return performAction(ACTION_ALL & ~ACTION_WRITE);
635 }
636
637 return true;
638}
639
640
642{
643 return performAction(ACTION_ALL);
644}
645
646
648{
649 if (&mpm.mesh() == &mesh_)
650 {
651 expire();
652 }
653
654 // pointMesh and interpolation will have been reset in mesh.update
655}
656
657
659{
660 if (&mesh == &mesh_)
661 {
662 expire();
663 }
664}
665
666
668{
669 if (state != polyMesh::UNCHANGED)
670 {
671 // May want to use force expiration here
672 expire();
673 }
674}
675
676
677bool Foam::sampledSurfaces::needsUpdate() const
678{
679 for (const sampledSurface& s : surfaces())
680 {
681 if (s.needsUpdate())
682 {
683 return true;
684 }
685 }
686
687 return false;
688}
689
690
691bool Foam::sampledSurfaces::expire(const bool force)
692{
693 // Dimension as fraction of mesh bounding box
694 const scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
695
696 label nChanged = 0;
697
698 forAll(*this, surfi)
699 {
700 sampledSurface& s = (*this)[surfi];
701
702 if (s.invariant() && !force)
703 {
704 // 'Invariant' - does not change when geometry does
705 continue;
706 }
707 if (s.expire())
708 {
709 ++nChanged;
710 }
711
712 writers_[surfi].expire();
713 writers_[surfi].mergeDim(mergeDim);
714 nFaces_[surfi] = 0;
715 }
716
717 // True if any surfaces just expired
718 return nChanged;
719}
720
721
722bool Foam::sampledSurfaces::update()
723{
724 if (!needsUpdate())
725 {
726 return false;
727 }
728
729 label nUpdated = 0;
730
731 forAll(*this, surfi)
732 {
733 sampledSurface& s = (*this)[surfi];
734
735 if (s.update())
736 {
737 ++nUpdated;
738 writers_[surfi].expire();
739 }
740
741 nFaces_[surfi] = returnReduce(s.faces().size(), sumOp<label>());
742 }
743
744 return nUpdated;
745}
746
747
749{
750 return mergeTol_;
751}
752
753
754Foam::scalar Foam::sampledSurfaces::mergeTol(const scalar tol) noexcept
755{
756 const scalar old(mergeTol_);
757 mergeTol_ = tol;
758 return old;
759}
760
761
762// ************************************************************************* //
Various functions to operate on Lists.
bool found
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.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static void mapCombineAllGather(const List< commsStruct > &comms, Container &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:97
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool merge(const dictionary &dict)
Merge entries from the given dictionary.
Definition: dictionary.C:812
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find for an entry (non-const access) with the given keyword.
Definition: dictionaryI.H:97
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
int verbose() const noexcept
Output verbosity level.
Definition: ensightMesh.C:125
bool expire()
Mark as needing an update.
Definition: ensightMeshI.H:36
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
static bool clean(std::string &str)
Definition: fileName.C:199
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:675
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:71
An abstract class for surfaces with sampling.
Set of surfaces to sample.
virtual bool read(const dictionary &dict)
Read the sampledSurfaces dictionary.
static scalar mergeTol() noexcept
Get merge tolerance.
virtual bool execute()
Sample and store if the sampleOnExecute is enabled.
virtual bool write()
Sample and write.
splitCell * master() const
Definition: splitCell.H:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
engineTime & runTime
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
bool found(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as found_if.
Definition: ListOps.H:679
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:55
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const direction noexcept
Definition: Scalar.H:223
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Foam::surfaceFields.
const dictionary formatOptions(propsDict.subOrEmptyDict("formatOptions", keyType::LITERAL))