modifyMesh.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-2021 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
27Application
28 modifyMesh
29
30Group
31 grpMeshAdvancedUtilities
32
33Description
34 Manipulate mesh elements.
35
36 Actions are:
37 (boundary)points:
38 - move
39
40 (boundary)edges:
41 - split and move introduced point
42
43 (boundary)faces:
44 - split(triangulate) and move introduced point
45
46 edges:
47 - collapse
48
49 cells:
50 - split into polygonal base pyramids around newly introduced mid
51 point
52
53 Is a bit of a loose collection of mesh change drivers.
54
55\*---------------------------------------------------------------------------*/
56
57#include "argList.H"
58#include "Time.H"
59#include "polyMesh.H"
60#include "polyTopoChange.H"
61#include "mapPolyMesh.H"
62#include "boundaryCutter.H"
63#include "cellSplitter.H"
64#include "edgeCollapser.H"
65#include "meshTools.H"
66#include "Pair.H"
67#include "globalIndex.H"
68#include "topoSet.H"
69#include "processorMeshes.H"
70#include "IOdictionary.H"
71
72using namespace Foam;
73
74// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75
76// Locate point on patch. Returns (mesh) point label.
77label findPoint(const primitivePatch& pp, const point& nearPoint)
78{
79 const pointField& points = pp.points();
80 const labelList& meshPoints = pp.meshPoints();
81
82 // Find nearest and next nearest
83 scalar minDistSqr = GREAT;
84 label minI = -1;
85
86 scalar almostMinDistSqr = GREAT;
87 label almostMinI = -1;
88
89 for (const label pointi : meshPoints)
90 {
91 scalar distSqr = magSqr(nearPoint - points[pointi]);
92
93 if (distSqr < minDistSqr)
94 {
95 almostMinDistSqr = minDistSqr;
96 almostMinI = minI;
97
98 minDistSqr = distSqr;
99 minI = pointi;
100 }
101 else if (distSqr < almostMinDistSqr)
102 {
103 almostMinDistSqr = distSqr;
104 almostMinI = pointi;
105 }
106 }
107
108
109 // Decide if nearPoint unique enough.
110 Info<< "Found to point " << nearPoint << nl
111 << " nearest point : " << minI
112 << " distance " << Foam::sqrt(minDistSqr)
113 << " at " << points[minI] << nl
114 << " next nearest point : " << almostMinI
115 << " distance " << Foam::sqrt(almostMinDistSqr)
116 << " at " << points[almostMinI] << endl;
117
118 if (almostMinDistSqr < 4*minDistSqr)
119 {
120 Info<< "Next nearest too close to nearest. Aborting" << endl;
121
122 return -1;
123 }
124 else
125 {
126 return minI;
127 }
128}
129
130
131// Locate edge on patch. Return mesh edge label.
132label findEdge
133(
134 const primitiveMesh& mesh,
135 const primitivePatch& pp,
136 const point& nearPoint
137)
138{
139 const pointField& localPoints = pp.localPoints();
140 const pointField& points = pp.points();
141 const labelList& meshPoints = pp.meshPoints();
142 const edgeList& edges = pp.edges();
143
144 // Find nearest and next nearest
145 scalar minDist = GREAT;
146 label minI = -1;
147
148 scalar almostMinDist = GREAT;
149 label almostMinI = -1;
150
151 for (const edge& e : edges)
152 {
153 pointHit pHit(e.line(localPoints).nearestDist(nearPoint));
154
155 if (pHit.hit())
156 {
157 if (pHit.distance() < minDist)
158 {
159 almostMinDist = minDist;
160 almostMinI = minI;
161
162 minDist = pHit.distance();
163 minI = meshTools::findEdge
164 (
165 mesh,
166 meshPoints[e[0]],
167 meshPoints[e[1]]
168 );
169 }
170 else if (pHit.distance() < almostMinDist)
171 {
172 almostMinDist = pHit.distance();
173 almostMinI = meshTools::findEdge
174 (
175 mesh,
176 meshPoints[e[0]],
177 meshPoints[e[1]]
178 );
179 }
180 }
181 }
182
183 if (minI == -1)
184 {
185 Info<< "Did not find edge close to point " << nearPoint << endl;
186
187 return -1;
188 }
189
190
191 // Decide if nearPoint unique enough.
192 Info<< "Found to point " << nearPoint << nl
193 << " nearest edge : " << minI
194 << " distance " << minDist << " endpoints "
195 << mesh.edges()[minI].line(points) << nl
196 << " next nearest edge : " << almostMinI
197 << " distance " << almostMinDist << " endpoints "
198 << mesh.edges()[almostMinI].line(points) << nl
199 << endl;
200
201 if (almostMinDist < 2*minDist)
202 {
203 Info<< "Next nearest too close to nearest. Aborting" << endl;
204
205 return -1;
206 }
207 else
208 {
209 return minI;
210 }
211}
212
213
214// Find face on patch. Return mesh face label.
215label findFace
216(
217 const primitiveMesh& mesh,
218 const primitivePatch& pp,
219 const point& nearPoint
220)
221{
222 const pointField& points = pp.points();
223
224 // Find nearest and next nearest
225 scalar minDist = GREAT;
226 label minI = -1;
227
228 scalar almostMinDist = GREAT;
229 label almostMinI = -1;
230
231 forAll(pp, patchFacei)
232 {
233 pointHit pHit(pp[patchFacei].nearestPoint(nearPoint, points));
234
235 if (pHit.hit())
236 {
237 if (pHit.distance() < minDist)
238 {
239 almostMinDist = minDist;
240 almostMinI = minI;
241
242 minDist = pHit.distance();
243 minI = patchFacei + mesh.nInternalFaces();
244 }
245 else if (pHit.distance() < almostMinDist)
246 {
247 almostMinDist = pHit.distance();
248 almostMinI = patchFacei + mesh.nInternalFaces();
249 }
250 }
251 }
252
253 if (minI == -1)
254 {
255 Info<< "Did not find face close to point " << nearPoint << endl;
256
257 return -1;
258 }
259
260
261 // Decide if nearPoint unique enough.
262 Info<< "Found to point " << nearPoint << nl
263 << " nearest face : " << minI
264 << " distance " << minDist
265 << " to face centre " << mesh.faceCentres()[minI] << nl
266 << " next nearest face : " << almostMinI
267 << " distance " << almostMinDist
268 << " to face centre " << mesh.faceCentres()[almostMinI] << nl
269 << endl;
270
271 if (almostMinDist < 2*minDist)
272 {
273 Info<< "Next nearest too close to nearest. Aborting" << endl;
274
275 return -1;
276 }
277 else
278 {
279 return minI;
280 }
281}
282
283
284// Find cell with cell centre close to given point.
285label findCell(const primitiveMesh& mesh, const point& nearPoint)
286{
287 label celli = mesh.findCell(nearPoint);
288
289 if (celli != -1)
290 {
291 scalar distToCcSqr = magSqr(nearPoint - mesh.cellCentres()[celli]);
292
293 const labelList& cPoints = mesh.cellPoints()[celli];
294
295 label minI = -1;
296 scalar minDistSqr = GREAT;
297
298 for (const label pointi : cPoints)
299 {
300 scalar distSqr = magSqr(nearPoint - mesh.points()[pointi]);
301
302 if (distSqr < minDistSqr)
303 {
304 minDistSqr = distSqr;
305 minI = pointi;
306 }
307 }
308
309 // Decide if nearPoint unique enough.
310 Info<< "Found to point " << nearPoint << nl
311 << " nearest cell : " << celli
312 << " distance " << Foam::sqrt(distToCcSqr)
313 << " to cell centre " << mesh.cellCentres()[celli] << nl
314 << " nearest mesh point : " << minI
315 << " distance " << Foam::sqrt(minDistSqr)
316 << " to " << mesh.points()[minI] << nl
317 << endl;
318
319 if (minDistSqr < 4*distToCcSqr)
320 {
321 Info<< "Mesh point too close to nearest cell centre. Aborting"
322 << endl;
323
324 celli = -1;
325 }
326 }
327
328 return celli;
329}
330
331
332
333
334int main(int argc, char *argv[])
335{
336 argList::addNote
337 (
338 "Manipulate mesh elements.\n"
339 "For example, moving points, splitting/collapsing edges etc."
340 );
341 #include "addOverwriteOption.H"
342 argList::addOption("dict", "file", "Alternative modifyMeshDict");
343
344 argList::noFunctionObjects(); // Never use function objects
345
346 #include "setRootCase.H"
347 #include "createTime.H"
348 #include "createPolyMesh.H"
349
350 const word oldInstance = mesh.pointsInstance();
351
352 const bool overwrite = args.found("overwrite");
353
354 Info<< "Reading modifyMeshDict\n" << endl;
355
356 // Read meshing dictionary
357 const word dictName("modifyMeshDict");
359 const IOdictionary dict(dictIO);
360
361 // Read all from the dictionary.
362 List<Pair<point>> pointsToMove(dict.lookup("pointsToMove"));
363 List<Pair<point>> edgesToSplit(dict.lookup("edgesToSplit"));
364 List<Pair<point>> facesToTriangulate
365 (
366 dict.lookup("facesToTriangulate")
367 );
368 // Optional
369 List<Pair<point>> facesToSplit;
370 dict.readIfPresent("facesToSplit", facesToSplit);
371
372 bool cutBoundary =
373 (
374 pointsToMove.size()
375 || edgesToSplit.size()
376 || facesToTriangulate.size()
377 || facesToSplit.size()
378 );
379
380 List<Pair<point>> edgesToCollapse(dict.lookup("edgesToCollapse"));
381
382 bool collapseEdge = edgesToCollapse.size();
383
384 List<Pair<point>> cellsToPyramidise(dict.lookup("cellsToSplit"));
385
386 bool cellsToSplit = cellsToPyramidise.size();
387
388 // List<Tuple2<pointField,point>>
389 // cellsToCreate(dict.lookup("cellsToCreate"));
390
391 Info<< "Read from " << dict.name() << nl
392 << " Boundary cutting module:" << nl
393 << " points to move :" << pointsToMove.size() << nl
394 << " edges to split :" << edgesToSplit.size() << nl
395 << " faces to split :" << facesToSplit.size() << nl
396 << " faces to triangulate:" << facesToTriangulate.size() << nl
397 << " Cell splitting module:" << nl
398 << " cells to split :" << cellsToPyramidise.size() << nl
399 << " Edge collapsing module:" << nl
400 << " edges to collapse :" << edgesToCollapse.size() << nl
401 //<< " cells to create :" << cellsToCreate.size() << nl
402 << endl;
403
404 if
405 (
406 (cutBoundary && collapseEdge)
407 || (cutBoundary && cellsToSplit)
408 || (collapseEdge && cellsToSplit)
409 )
410 {
412 << "Used more than one mesh modifying module "
413 << "(boundary cutting, cell splitting, edge collapsing)" << nl
414 << "Please do them in separate passes." << exit(FatalError);
415 }
416
417
418
419 // Get calculating engine for all of outside
420 const SubList<face> outsideFaces
421 (
422 mesh.faces(),
423 mesh.nBoundaryFaces(),
424 mesh.nInternalFaces()
425 );
426
427 primitivePatch allBoundary(outsideFaces, mesh.points());
428
429
430 // Look up mesh labels and convert to input for boundaryCutter.
431
432 bool validInputs = true;
433
434
435 Info<< nl << "Looking up points to move ..." << nl << endl;
436 Map<point> pointToPos(pointsToMove.size());
437 for (const Pair<point>& pts : pointsToMove)
438 {
439 const label pointi = findPoint(allBoundary, pts.first());
440
441 if (pointi == -1 || !pointToPos.insert(pointi, pts.second()))
442 {
443 Info<< "Could not insert mesh point " << pointi
444 << " for input point " << pts.first() << nl
445 << "Perhaps the point is already marked for moving?" << endl;
446 validInputs = false;
447 }
448 }
449
450
451 Info<< nl << "Looking up edges to split ..." << nl << endl;
452 Map<List<point>> edgeToCuts(edgesToSplit.size());
453 for (const Pair<point>& pts : edgesToSplit)
454 {
455 label edgeI = findEdge(mesh, allBoundary, pts.first());
456
457 if
458 (
459 edgeI == -1
460 || !edgeToCuts.insert(edgeI, List<point>(1, pts.second()))
461 )
462 {
463 Info<< "Could not insert mesh edge " << edgeI
464 << " for input point " << pts.first() << nl
465 << "Perhaps the edge is already marked for cutting?" << endl;
466
467 validInputs = false;
468 }
469 }
470
471
472 Info<< nl << "Looking up faces to triangulate ..." << nl << endl;
473 Map<point> faceToDecompose(facesToTriangulate.size());
474 for (const Pair<point>& pts : facesToTriangulate)
475 {
476 label facei = findFace(mesh, allBoundary, pts.first());
477
478 if (facei == -1 || !faceToDecompose.insert(facei, pts.second()))
479 {
480 Info<< "Could not insert mesh face " << facei
481 << " for input point " << pts.first() << nl
482 << "Perhaps the face is already marked for splitting?" << endl;
483
484 validInputs = false;
485 }
486 }
487
488
489 Info<< nl << "Looking up faces to split ..." << nl << endl;
490 Map<labelPair> faceToSplit(facesToSplit.size());
491 for (const Pair<point>& pts : facesToSplit)
492 {
493 label facei = findFace(mesh, allBoundary, pts.first());
494 if (facei == -1)
495 {
496 Info<< "Could not insert mesh face " << facei
497 << " for input point " << pts.first() << nl
498 << "Perhaps the face is already marked for splitting?" << endl;
499
500 validInputs = false;
501 }
502 else
503 {
504 // Find nearest points on face
505 const primitivePatch pp
506 (
507 SubList<face>(mesh.faces(), 1, facei),
508 mesh.points()
509 );
510
511 const label p0 = findPoint(pp, pts.first());
512 const label p1 = findPoint(pp, pts.second());
513
514 const face& f = mesh.faces()[facei];
515
516 if
517 (
518 p0 != -1
519 && p1 != -1
520 && faceToSplit.insert(facei, labelPair(f.find(p0), f.find(p1)))
521 )
522 {}
523 else
524 {
525 Info<< "Could not insert mesh face " << facei
526 << " for input coordinates " << pts
527 << " with vertices " << p0 << " and " << p1 << nl
528 << "Perhaps the face is already marked for splitting?"
529 << endl;
530
531 }
532 }
533 }
534
535
536 Info<< nl << "Looking up cells to convert to pyramids around"
537 << " cell centre ..." << nl << endl;
538 Map<point> cellToPyrCentre(cellsToPyramidise.size());
539 for (const Pair<point>& pts : cellsToPyramidise)
540 {
541 label celli = findCell(mesh, pts.first());
542
543 if (celli == -1 || !cellToPyrCentre.insert(celli, pts.second()))
544 {
545 Info<< "Could not insert mesh cell " << celli
546 << " for input point " << pts.first() << nl
547 << "Perhaps the cell is already marked for splitting?" << endl;
548
549 validInputs = false;
550 }
551 }
552
553
554 Info<< nl << "Looking up edges to collapse ..." << nl << endl;
555 Map<point> edgeToPos(edgesToCollapse.size());
556 for (const Pair<point>& pts : edgesToCollapse)
557 {
558 label edgeI = findEdge(mesh, allBoundary, pts.first());
559
560 if (edgeI == -1 || !edgeToPos.insert(edgeI, pts.second()))
561 {
562 Info<< "Could not insert mesh edge " << edgeI
563 << " for input point " << pts.first() << nl
564 << "Perhaps the edge is already marked for collaping?" << endl;
565
566 validInputs = false;
567 }
568 }
569
570
571
572 if (!validInputs)
573 {
574 Info<< nl << "There was a problem in one of the inputs in the"
575 << " dictionary. Not modifying mesh." << endl;
576 }
577 else if (cellToPyrCentre.size())
578 {
579 Info<< nl << "All input cells located. Modifying mesh." << endl;
580
581 // Mesh change engine
582 cellSplitter cutter(mesh);
583
584 // Topo change container
585 polyTopoChange meshMod(mesh);
586
587 // Insert commands into meshMod
588 cutter.setRefinement(cellToPyrCentre, meshMod);
589
590 // Do changes
591 autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
592
593 if (morphMap().hasMotionPoints())
594 {
595 mesh.movePoints(morphMap().preMotionPoints());
596 }
597
598 cutter.updateMesh(morphMap());
599
600 if (!overwrite)
601 {
602 ++runTime;
603 }
604 else
605 {
606 mesh.setInstance(oldInstance);
607 }
608
609 // Write resulting mesh
610 Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
611 mesh.write();
612 topoSet::removeFiles(mesh);
613 processorMeshes::removeFiles(mesh);
614 }
615 else if (edgeToPos.size())
616 {
617 Info<< nl << "All input edges located. Modifying mesh." << endl;
618
619 // Mesh change engine
620 edgeCollapser cutter(mesh);
621
622 const edgeList& edges = mesh.edges();
623 const pointField& points = mesh.points();
624
625 pointField newPoints(points);
626
627 bitSet collapseEdge(mesh.nEdges());
628 Map<point> collapsePointToLocation(mesh.nPoints());
629
630 // Get new positions and construct collapse network
631 forAllConstIters(edgeToPos, iter)
632 {
633 label edgeI = iter.key();
634 const edge& e = edges[edgeI];
635
636 collapseEdge.set(edgeI);
637 collapsePointToLocation.set(e[1], points[e[0]]);
638
639 newPoints[e[0]] = iter.val();
640 }
641
642 // Move master point to destination.
643 mesh.movePoints(newPoints);
644
645 List<pointEdgeCollapse> allPointInfo;
646 const globalIndex globalPoints(mesh.nPoints());
647 labelList pointPriority(mesh.nPoints(), Zero);
648
649 cutter.consistentCollapse
650 (
652 pointPriority,
653 collapsePointToLocation,
655 allPointInfo
656 );
657
658 // Topo change container
659 polyTopoChange meshMod(mesh);
660
661 // Insert
662 cutter.setRefinement(allPointInfo, meshMod);
663
664 // Do changes
665 autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
666
667 if (morphMap().hasMotionPoints())
668 {
669 mesh.movePoints(morphMap().preMotionPoints());
670 }
671
672 // Not implemented yet:
673 //cutter.updateMesh(morphMap());
674
675
676 if (!overwrite)
677 {
678 ++runTime;
679 }
680 else
681 {
682 mesh.setInstance(oldInstance);
683 }
684
685 // Write resulting mesh
686 Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
687 mesh.write();
688 topoSet::removeFiles(mesh);
689 processorMeshes::removeFiles(mesh);
690 }
691 else
692 {
693 Info<< nl << "All input points located. Modifying mesh." << endl;
694
695 // Mesh change engine
696 boundaryCutter cutter(mesh);
697
698 // Topo change container
699 polyTopoChange meshMod(mesh);
700
701 // Insert commands into meshMod
702 cutter.setRefinement
703 (
704 pointToPos,
705 edgeToCuts,
706 faceToSplit, // Faces to split diagonally
707 faceToDecompose, // Faces to triangulate
708 meshMod
709 );
710
711 // Do changes
712 autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
713
714 if (morphMap().hasMotionPoints())
715 {
716 mesh.movePoints(morphMap().preMotionPoints());
717 }
718
719 cutter.updateMesh(morphMap());
720
721 if (!overwrite)
722 {
723 ++runTime;
724 }
725 else
726 {
727 mesh.setInstance(oldInstance);
728 }
729
730 // Write resulting mesh
731 Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
732 mesh.write();
733 topoSet::removeFiles(mesh);
734 processorMeshes::removeFiles(mesh);
735 }
736
737
738 Info<< "\nEnd\n" << endl;
739 return 0;
740}
741
742
743// ************************************************************************* //
static constexpr label size() noexcept
Return the number of elements in the FixedList.
Definition: FixedList.H:416
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
A list of faces which address into the list of points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
A List obtained as a section of another List.
Definition: SubList.H:70
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
Does modifications to boundary faces.
Does pyramidal decomposition of selected cells. So all faces will become base of pyramid with as top ...
Definition: cellSplitter.H:60
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:69
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:103
Direct mesh changes based on v1.3 polyTopoChange syntax.
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
A class for handling words, derived from Foam::string.
Definition: word.H:68
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
const volScalarField & p0
Definition: EEqn.H:36
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const word dictName("faMeshDefinition")
const pointField & points
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:359
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
IOobject dictIO
Foam::argList args(argc, argv)
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