interfaceTrackingFvMesh.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) 2019 Zeljko Tukovic, FSB Zagreb.
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
31 #include "motionSolver.H"
32 #include "volFields.H"
33 #include "wedgeFaPatch.H"
34 #include "wedgeFaPatchFields.H"
35 #include "slipFaPatchFields.H"
37 #include "slipFvPatchFields.H"
38 #include "symmetryFvPatchFields.H"
39 #include "wallFvPatch.H"
40 #include "polyPatchID.H"
41 #include "fvcMeshPhi.H"
43 #include "EulerDdtScheme.H"
44 #include "CrankNicolsonDdtScheme.H"
45 #include "backwardDdtScheme.H"
46 #include "twoDPointCorrector.H"
47 #include "gravityMeshObject.H"
49 #include "demandDrivenData.H"
50 #include "unitConversion.H"
51 
52 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
58  (
61  IOobject
62  );
63 }
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
68 void Foam::interfaceTrackingFvMesh::initializeData()
69 {
70  // Set free surface patch index
71  {
72  const word fsPatchName(motion().get<word>("fsPatchName"));
73 
74  polyPatchID patch(fsPatchName, this->boundaryMesh());
75 
76  if (!patch.active())
77  {
79  << "Patch name " << fsPatchName << " not found."
80  << abort(FatalError);
81  }
82 
83  fsPatchIndex_ = patch.index();
84  }
85 
86  // Set point normal correction for finite area mesh
87  {
89 
90  for (const word& patchName : pointNormalsCorrectionPatches_)
91  {
92  label patchID = aMesh().boundary().findPatchID(patchName);
93 
94  if (patchID == -1)
95  {
97  << "Patch name '" << patchName
98  << "' for point normals correction does not exist"
99  << abort(FatalError);
100  }
101 
102  correction[patchID] = true;
103  }
104  }
105 
106  // Read motion direction
107  if (!normalMotionDir_)
108  {
109  motionDir_ = normalised(motion().get<vector>("motionDir"));
110  }
111 
112  // Check if contact angle is defined
113  makeContactAngle();
114 
116  (
117  "nonReflectingFreeSurfacePatches",
118  nonReflectingFreeSurfacePatches_
119  );
120 }
121 
122 
123 void Foam::interfaceTrackingFvMesh::makeUs() const
124 {
126  << "making free-surface velocity field" << nl;
127 
128  if (UsPtr_)
129  {
131  << "free-surface velocity field already exists"
132  << abort(FatalError);
133  }
134 
135  wordList patchFieldTypes
136  (
137  aMesh().boundary().size(),
138  zeroGradientFaPatchVectorField::typeName
139  );
140 
141  forAll(aMesh().boundary(), patchI)
142  {
143  if
144  (
145  aMesh().boundary()[patchI].type()
146  == wedgeFaPatch::typeName
147  )
148  {
149  patchFieldTypes[patchI] =
150  wedgeFaPatchVectorField::typeName;
151  }
152  else
153  {
154  label ngbPolyPatchID =
155  aMesh().boundary()[patchI].ngbPolyPatchIndex();
156 
157  if (ngbPolyPatchID != -1)
158  {
159  if
160  (
161  mesh().boundary()[ngbPolyPatchID].type()
162  == wallFvPatch::typeName
163  )
164  {
165  patchFieldTypes[patchI] =
166  slipFaPatchVectorField::typeName;
167  }
168  }
169  }
170  }
171 
172  for (const word& patchName : fixedFreeSurfacePatches_)
173  {
174  const label fixedPatchID =
175  aMesh().boundary().findPatchID(patchName);
176 
177  if (fixedPatchID == -1)
178  {
180  << "Wrong faPatch name '" << patchName
181  << "' in the fixedFreeSurfacePatches list"
182  << " defined in the dynamicMeshDict dictionary"
183  << abort(FatalError);
184  }
185 
186  label ngbPolyPatchID =
187  aMesh().boundary()[fixedPatchID].ngbPolyPatchIndex();
188 
189  if (ngbPolyPatchID != -1)
190  {
191  if
192  (
193  mesh().boundary()[ngbPolyPatchID].type()
194  == wallFvPatch::typeName
195  )
196  {
197  patchFieldTypes[fixedPatchID] =
198  fixedValueFaPatchVectorField::typeName;
199  }
200  }
201  }
202 
203  UsPtr_ = new areaVectorField
204  (
205  IOobject
206  (
207  "Us",
208  mesh().time().timeName(),
209  mesh(),
212  ),
213  aMesh(),
215  patchFieldTypes
216  );
217 
218  for (const word& patchName : fixedFreeSurfacePatches_)
219  {
220  const label fixedPatchID = aMesh().boundary().findPatchID(patchName);
221 
222  if (fixedPatchID == -1)
223  {
225  << "Wrong faPatch name '" << patchName
226  << "' in the fixedFreeSurfacePatches list"
227  << " defined in the dynamicMeshDict dictionary"
228  << abort(FatalError);
229  }
230 
231  label ngbPolyPatchID =
232  aMesh().boundary()[fixedPatchID].ngbPolyPatchIndex();
233 
234  if (ngbPolyPatchID != -1)
235  {
236  if
237  (
238  mesh().boundary()[ngbPolyPatchID].type()
239  == wallFvPatch::typeName
240  )
241  {
242  UsPtr_->boundaryFieldRef()[fixedPatchID] == Zero;
243  }
244  }
245  }
246 }
247 
248 
249 void Foam::interfaceTrackingFvMesh::makeFsNetPhi() const
250 {
252  << "making free-surface net flux" << nl;
253 
254  if (fsNetPhiPtr_)
255  {
257  << "free-surface net flux already exists"
258  << abort(FatalError);
259  }
260 
261  fsNetPhiPtr_ = new areaScalarField
262  (
263  IOobject
264  (
265  "fsNetPhi",
266  mesh().time().timeName(),
267  mesh(),
270  ),
271  aMesh(),
273  );
274 }
275 
276 
277 void Foam::interfaceTrackingFvMesh::makeControlPoints()
278 {
280  << "making control points" << nl;
281 
282  if (controlPointsPtr_)
283  {
285  << "control points already exists"
286  << abort(FatalError);
287  }
288 
289  IOobject controlPointsHeader
290  (
291  "controlPoints",
292  mesh().time().timeName(),
293  mesh(),
295  );
296 
297  if (controlPointsHeader.typeHeaderOk<vectorIOField>())
298  {
299  Info<< "Reading control points" << endl;
300  controlPointsPtr_ =
301  new vectorIOField
302  (
303  IOobject
304  (
305  "controlPoints",
306  mesh().time().timeName(),
307  mesh(),
310  )
311  );
312  }
313  else
314  {
315  Info<< "Creating new control points" << endl;
316  controlPointsPtr_ =
317  new vectorIOField
318  (
319  IOobject
320  (
321  "controlPoints",
322  mesh().time().timeName(),
323  mesh(),
326  ),
327  aMesh().areaCentres().internalField()
328  );
329 
330  initializeControlPointsPosition();
331  }
332 }
333 
334 
335 void Foam::interfaceTrackingFvMesh::makeMotionPointsMask()
336 {
338  << "making motion points mask" << nl;
339 
340  if (motionPointsMaskPtr_)
341  {
343  << "motion points mask already exists"
344  << abort(FatalError);
345  }
346 
347  motionPointsMaskPtr_ = new labelList
348  (
349  mesh().boundaryMesh()[fsPatchIndex()].nPoints(),
350  1
351  );
352 
353  // Mark free surface boundary points
354  // that belong to processor patches
355  forAll(aMesh().boundary(), patchI)
356  {
357  if
358  (
359  aMesh().boundary()[patchI].type()
360  == processorFaPatch::typeName
361  )
362  {
363  const labelList& patchPoints =
364  aMesh().boundary()[patchI].pointLabels();
365 
366  forAll(patchPoints, pointI)
367  {
368  motionPointsMask()[patchPoints[pointI]] = -1;
369  }
370  }
371  }
372 
373  // Mark fixed free surface boundary points
374  for (const word& patchName : fixedFreeSurfacePatches_)
375  {
376  const label fixedPatchID = aMesh().boundary().findPatchID(patchName);
377 
378  if (fixedPatchID == -1)
379  {
381  << "Wrong faPatch name in the fixedFreeSurfacePatches list"
382  << " defined in the dynamicMeshDict dictionary"
383  << abort(FatalError);
384  }
385 
386  const labelList& patchPoints =
387  aMesh().boundary()[fixedPatchID].pointLabels();
388 
389  forAll(patchPoints, pointI)
390  {
391  motionPointsMask()[patchPoints[pointI]] = 0;
392  }
393  }
394 }
395 
396 
397 void Foam::interfaceTrackingFvMesh::makeDirections()
398 {
400  << "make displacement directions for points and control points" << nl;
401 
402  if (pointsDisplacementDirPtr_ || facesDisplacementDirPtr_)
403  {
405  << "points, control points displacement directions already exist"
406  << abort(FatalError);
407  }
408 
409  pointsDisplacementDirPtr_ =
410  new vectorField
411  (
412  mesh().boundaryMesh()[fsPatchIndex()].nPoints(),
413  Zero
414  );
415 
416  facesDisplacementDirPtr_ =
417  new vectorField
418  (
419  mesh().boundaryMesh()[fsPatchIndex()].size(),
420  Zero
421  );
422 
423  if (!normalMotionDir())
424  {
425  if (mag(motionDir_) < SMALL)
426  {
428  << "Zero motion direction"
429  << abort(FatalError);
430  }
431 
432  facesDisplacementDir() = motionDir_;
433  pointsDisplacementDir() = motionDir_;
434  }
435 
436  updateDisplacementDirections();
437 }
438 
439 
440 void Foam::interfaceTrackingFvMesh::makePhis()
441 {
443  << "making free-surface flux" << nl;
444 
445  if (phisPtr_)
446  {
448  << "free-surface flux already exists"
449  << abort(FatalError);
450  }
451 
452  phisPtr_ = new edgeScalarField
453  (
454  IOobject
455  (
456  "phis",
457  mesh().time().timeName(),
458  mesh(),
461  ),
462  linearEdgeInterpolate(Us()) & aMesh().Le()
463  );
464 }
465 
466 
467 void Foam::interfaceTrackingFvMesh::makeSurfactConc() const
468 {
470  << "making free-surface surfactant concentration field" << nl;
471 
472  if (surfactConcPtr_)
473  {
475  << "free-surface surfactant concentration field already exists"
476  << abort(FatalError);
477  }
478 
479  surfactConcPtr_ = new areaScalarField
480  (
481  IOobject
482  (
483  "Cs",
484  mesh().time().timeName
485  (
486  mesh().time().startTime().value()
487  ),
488  // mesh().time().timeName(),
489  mesh(),
492  ),
493  aMesh()
494  );
495 }
496 
497 
498 void Foam::interfaceTrackingFvMesh::makeBulkSurfactConc() const
499 {
501  << "making volume surfactant concentration field" << nl;
502 
503  if (bulkSurfactConcPtr_)
504  {
506  << "volume surfactant concentration field already exists"
507  << abort(FatalError);
508  }
509 
510  bulkSurfactConcPtr_ = new volScalarField
511  (
512  IOobject
513  (
514  "C",
515  mesh().time().timeName
516  (
517  mesh().time().startTime().value()
518  ),
519  // mesh().time().timeName(),
520  mesh(),
523  ),
524  mesh()
525  );
526  volScalarField& bulkSurfactConc = *bulkSurfactConcPtr_;
527 
528  if (mesh().time().timeIndex()-1 == 0)
529  {
530  // Initialize uniform volume surfactant concentration
531  bulkSurfactConc = surfactant().bulkConc();
532  bulkSurfactConc.correctBoundaryConditions();
533  }
534 }
535 
536 
537 void Foam::interfaceTrackingFvMesh::makeSurfaceTension() const
538 {
540  << "making surface tension field" << nl;
541 
542  if (surfaceTensionPtr_)
543  {
545  << "surface tension field already exists"
546  << abort(FatalError);
547  }
548 
549  surfaceTensionPtr_ = new areaScalarField
550  (
551  IOobject
552  (
553  "surfaceTension",
554  mesh().time().timeName(),
555  mesh(),
558  ),
559  sigma() + surfactant().dSigma(surfactantConcentration())/rho_
560  );
561 }
562 
563 
564 void Foam::interfaceTrackingFvMesh::makeSurfactant() const
565 {
567  << "making surfactant properties" << nl;
568 
569  if (surfactantPtr_)
570  {
572  << "surfactant properties already exists"
573  << abort(FatalError);
574  }
575 
576  const dictionary& surfactProp =
577  motion().subDict("surfactantProperties");
578 
579  surfactantPtr_ = new surfactantProperties(surfactProp);
580 }
581 
582 
583 void Foam::interfaceTrackingFvMesh::makeContactAngle()
584 {
586  << "making contact angle field" << nl;
587 
588  if (contactAnglePtr_)
589  {
591  << "contact angle already exists"
592  << abort(FatalError);
593  }
594 
595  // Check if contactAngle is defined
596  IOobject contactAngleHeader
597  (
598  "contactAngle",
599  mesh().time().timeName(),
600  mesh(),
602  );
603 
604  if (contactAngleHeader.typeHeaderOk<areaScalarField>())
605  {
606  Info<< "Reading contact angle field" << endl;
607 
608  contactAnglePtr_ =
609  new areaScalarField
610  (
611  IOobject
612  (
613  "contactAngle",
614  mesh().time().timeName(),
615  mesh(),
618  ),
619  aMesh()
620  );
621  }
622 }
623 
624 
625 void Foam::interfaceTrackingFvMesh::updateDisplacementDirections()
626 {
627  if (normalMotionDir())
628  {
629  // Update point displacement direction
630  pointsDisplacementDir() = aMesh().pointAreaNormals();
631 
632  // Correct point displacement direction at contact line
633  forAll(aMesh().boundary(), patchI)
634  {
635  if (contactAnglePtr_)
636  {
637  label ngbPolyPatchID =
638  aMesh().boundary()[patchI].ngbPolyPatchIndex();
639 
640  if (ngbPolyPatchID != -1)
641  {
642  if
643  (
644  mesh().boundary()[ngbPolyPatchID].type()
645  == wallFvPatch::typeName
646  )
647  {
648  labelList patchPoints =
649  aMesh().boundary()[patchI].pointLabels();
650 
651  vectorField N
652  (
653  aMesh().boundary()[patchI]
654  .ngbPolyPatchPointNormals()
655  );
656 
657  forAll(patchPoints, pointI)
658  {
659  pointsDisplacementDir()[patchPoints[pointI]] -=
660  N[pointI]
661  *(
662  N[pointI]
663  & pointsDisplacementDir()[patchPoints[pointI]]
664  );
665 
666  pointsDisplacementDir()[patchPoints[pointI]] /=
667  mag
668  (
669  pointsDisplacementDir()
670  [
671  patchPoints[pointI]
672  ]
673  ) + SMALL;
674  }
675  }
676  }
677  }
678  }
679 
680  // Update face displacement direction
681  facesDisplacementDir() =
683 
684  // Correction of control points position
685  const vectorField& Cf = aMesh().areaCentres().internalField();
686 
687  controlPoints() =
688  facesDisplacementDir()
689  *(facesDisplacementDir()&(controlPoints() - Cf))
690  + Cf;
691  }
692 }
693 
694 
695 void Foam::interfaceTrackingFvMesh::initializeControlPointsPosition()
696 {
697  {
698  const faceList& faces = aMesh().faces();
699  const pointField& points = aMesh().points();
700 
701  pointField displacement(pointDisplacement());
702  scalarField sweptVolCorr(faces.size(), Zero);
703  correctPointDisplacement(sweptVolCorr, displacement);
704 
705  pointField newPoints(points + displacement);
706 
707  scalarField sweptVol(faces.size(), Zero);
708 
709  forAll(faces, faceI)
710  {
711  sweptVol[faceI] = -faces[faceI].sweptVol(points, newPoints);
712  }
713 
714  vectorField faceArea(faces.size(), Zero);
715 
716  forAll(faceArea, faceI)
717  {
718  faceArea[faceI] = faces[faceI].unitNormal(newPoints);
719  }
720 
721  scalarField deltaH = scalarField(aMesh().nFaces(), Zero);
722 
723  forAll(deltaH, faceI)
724  {
725  deltaH[faceI] = sweptVol[faceI]/
726  ((faceArea[faceI] & facesDisplacementDir()[faceI]) + SMALL);
727 
728  if (mag(faceArea[faceI] & facesDisplacementDir()[faceI]) < SMALL)
729  {
730  // Info<< (faceArea[faceI] & facesDisplacementDir()[faceI])
731  // << ", " << faceArea[faceI]
732  // << ", " << facesDisplacementDir()[faceI] << endl;
733 
734  FatalError
735  << "Something is wrong with specified motion direction"
736  << abort(FatalError);
737  }
738  }
739 
740  for (const word& patchName : fixedFreeSurfacePatches_)
741  {
742  label fixedPatchID = aMesh().boundary().findPatchID(patchName);
743 
744  if (fixedPatchID == -1)
745  {
746  FatalError
747  << "Wrong faPatch name in the fixedFreeSurfacePatches list"
748  << " defined in the freeSurfaceProperties dictionary"
749  << abort(FatalError);
750  }
751 
752  const labelList& eFaces =
753  aMesh().boundary()[fixedPatchID].edgeFaces();
754 
755  forAll(eFaces, edgeI)
756  {
757  deltaH[eFaces[edgeI]] *= 2.0;
758  }
759  }
760 
761  controlPoints() += facesDisplacementDir()*deltaH;
762  }
763 }
764 
765 
766 void Foam::interfaceTrackingFvMesh::smoothFreeSurfaceMesh()
767 {
768  Info<< "Smoothing free surface mesh" << endl;
769 
770  controlPoints() = aMesh().areaCentres().internalField();
771 
772  pointField displacement(pointDisplacement());
773 
774  const faceList& faces = aMesh().faces();
775  const pointField& points = aMesh().points();
776 
777  pointField newPoints(points + displacement);
778 
779  scalarField sweptVol(faces.size(), Zero);
780  forAll(faces, faceI)
781  {
782  sweptVol[faceI] = -faces[faceI].sweptVol(points, newPoints);
783  }
784 
785  vectorField faceArea(faces.size(), Zero);
786  forAll(faceArea, faceI)
787  {
788  faceArea[faceI] = faces[faceI].unitNormal(newPoints);
789  }
790 
791  scalarField deltaHf
792  (
793  sweptVol/(faceArea & facesDisplacementDir())
794  );
795 
796  for (const word& patchName : fixedFreeSurfacePatches_)
797  {
798  label fixedPatchID = aMesh().boundary().findPatchID(patchName);
799 
800  if (fixedPatchID == -1)
801  {
802  FatalError
803  << "Wrong faPatch name fixedFreeSurfacePatches list"
804  << " defined in the dynamicMeshDict dictionary"
805  << abort(FatalError);
806  }
807 
808  const labelList& eFaces =
809  aMesh().boundary()[fixedPatchID].edgeFaces();
810 
811  forAll(eFaces, edgeI)
812  {
813  deltaHf[eFaces[edgeI]] *= 2.0;
814  }
815  }
816 
817  controlPoints() += facesDisplacementDir()*deltaHf;
818 
819  displacement = pointDisplacement();
820 
821  velocityMotionSolver& vMotion =
822  refCast<velocityMotionSolver>
823  (
824  const_cast<motionSolver&>(motion())
825  );
826 
827  pointVectorField& pointMotionU = vMotion.pointMotionU();
828  pointMotionU.primitiveFieldRef() = Zero;
829 
830  fixedValuePointPatchVectorField& fsPatchPointMeshU =
831  refCast<fixedValuePointPatchVectorField>
832  (
833  const_cast<pointPatchVectorField&>
834  (
835  pointMotionU.boundaryField()[fsPatchIndex()]
836  )
837  );
838 
839  fsPatchPointMeshU ==
840  displacement/mesh().time().deltaT().value();
841 
843 }
844 
845 
846 void Foam::interfaceTrackingFvMesh::updateSurfaceFlux()
847 {
848  Phis() = fac::interpolate(Us()) & aMesh().Le();
849 }
850 
851 
852 void Foam::interfaceTrackingFvMesh::updateSurfactantConcentration()
853 {
854  if (!pureFreeSurface())
855  {
856  Info<< "Correct surfactant concentration" << endl << flush;
857 
858  updateSurfaceFlux();
859 
860  // Crate and solve the surfactanta transport equation
861  faScalarMatrix CsEqn
862  (
863  fam::ddt(surfactantConcentration())
864  + fam::div(Phis(), surfactantConcentration())
866  (
867  surfactant().diffusion(),
868  surfactantConcentration()
869  )
870  );
871 
872  if (surfactant().soluble())
873  {
874  #include "solveBulkSurfactant.H"
875 
877  (
878  IOobject
879  (
880  "Cb",
881  mesh().time().timeName(),
882  mesh(),
885  ),
886  aMesh(),
888  zeroGradientFaPatchScalarField::typeName
889  );
890 
891  Cb.ref().field() =
892  bulkSurfactantConcentration().boundaryField()[fsPatchIndex()];
893  Cb.correctBoundaryConditions();
894 
895  CsEqn +=
896  fam::Sp
897  (
898  surfactant().adsorptionCoeff()*Cb
899  + surfactant().adsorptionCoeff()
900  *surfactant().desorptionCoeff(),
901  surfactantConcentration()
902  )
903  - surfactant().adsorptionCoeff()
904  *Cb*surfactant().saturatedConc();
905  }
906 
907  CsEqn.solve();
908 
909  // Info<< "Correct surface tension" << endl;
910 
911  surfaceTension() =
912  sigma() + surfactant().dSigma(surfactantConcentration())/rho_;
913 
914  if (neg(min(surfaceTension().internalField().field())))
915  {
917  << "Surface tension is negative"
918  << abort(FatalError);
919  }
920  }
921 }
922 
923 
924 Foam::vector Foam::interfaceTrackingFvMesh::totalPressureForce() const
925 {
926  const scalarField& S = aMesh().S();
927 
929 
930  const scalarField& P = p().boundaryField()[fsPatchIndex()];
931 
932  vectorField pressureForces(S*P*n);
933 
934  return gSum(pressureForces);
935 }
936 
937 
938 Foam::vector Foam::interfaceTrackingFvMesh::totalViscousForce() const
939 {
940  const auto& turbulence =
941  mesh().lookupObject<turbulenceModel>("turbulenceProperties");
942 
943  scalarField nu(turbulence.nuEff(fsPatchIndex()));
944 
945  // const singlePhaseTransportModel& properties =
946  // mesh().lookupObject<singlePhaseTransportModel>
947  // (
948  // "transportProperties"
949  // );
950 
951  // dimensionedScalar nu("nu", properties);
952 
953  const scalarField& S = aMesh().S();
955 
956  vectorField nGradU
957  (
958  U().boundaryField()[fsPatchIndex()].snGrad()
959  );
960 
961  vectorField viscousForces
962  (
963  - nu*S
964  *(
965  nGradU
966  + (fac::grad(Us())().internalField()&n)
967  - (n*fac::div(Us())().internalField())
968  )
969  );
970 
971  return gSum(viscousForces);
972 }
973 
974 
975 Foam::vector Foam::interfaceTrackingFvMesh::totalSurfaceTensionForce() const
976 {
977  const scalarField& S = aMesh().S();
978 
980 
982 
983  vectorField surfTensionForces(n.size(), Zero);
984 
985  if (pureFreeSurface())
986  {
987  surfTensionForces =
988  S*sigma().value()
990  (
991  aMesh().Le()*aMesh().edgeLengthCorrection()
992  )().internalField();
993  }
994  else
995  {
996  surfTensionForces = surfaceTension().internalField().field()*K*S*n;
997  }
998 
999  return gSum(surfTensionForces);
1000 }
1001 
1002 
1003 Foam::scalar Foam::interfaceTrackingFvMesh::maxCourantNumber()
1004 {
1005  scalar CoNum = 0;
1006 
1007  if (pureFreeSurface())
1008  {
1009  const scalarField& dE = aMesh().lPN();
1010 
1011  CoNum = gMax
1012  (
1013  mesh().time().deltaT().value()/
1014  sqrt
1015  (
1016  Foam::pow(dE, 3.0)/2.0/M_PI/(sigma().value() + SMALL)
1017  )
1018  );
1019  }
1020  else
1021  {
1022  scalarField sigmaE
1023  (
1024  linearEdgeInterpolate(surfaceTension())().internalField().field()
1025  + SMALL
1026  );
1027 
1028  const scalarField& dE = aMesh().lPN();
1029 
1030  CoNum = gMax
1031  (
1032  mesh().time().deltaT().value()/
1033  sqrt
1034  (
1035  Foam::pow(dE, 3.0)/2.0/M_PI/sigmaE
1036  )
1037  );
1038  }
1039 
1040  return CoNum;
1041 }
1042 
1043 
1044 void Foam::interfaceTrackingFvMesh::updateProperties()
1045 {
1046  const singlePhaseTransportModel& properties =
1048  (
1049  "transportProperties"
1050  );
1051 
1052  rho_ = dimensionedScalar("rho", properties);
1053 
1054  sigma0_ = dimensionedScalar("sigma", properties)/rho_;
1055 }
1056 
1057 
1058 void Foam::interfaceTrackingFvMesh::correctPointDisplacement
1059 (
1060  const scalarField& sweptVolCorr,
1061  vectorField& displacement
1062 )
1063 {
1064  const labelListList& pFaces =
1065  aMesh().patch().pointFaces();
1066 
1067  const faceList& faces =
1068  aMesh().patch().localFaces();
1069 
1070  const pointField& points =
1071  aMesh().patch().localPoints();
1072 
1073  for (const word& patchName : fixedFreeSurfacePatches_)
1074  {
1075  label fixedPatchID = aMesh().boundary().findPatchID(patchName);
1076 
1077  const labelList& pLabels =
1078  aMesh().boundary()[fixedPatchID].pointLabels();
1079 
1080  const labelList& eFaces =
1081  aMesh().boundary()[fixedPatchID].edgeFaces();
1082 
1084 
1085  forAll(eFaces, edgeI)
1086  {
1087  label curFace = eFaces[edgeI];
1088 
1089  const labelList& curPoints = faces[curFace];
1090 
1091  forAll(curPoints, pointI)
1092  {
1093  label curPoint = curPoints[pointI];
1094  label index = pLabels.find(curPoint);
1095 
1096  if (index == -1)
1097  {
1098  if (!pointSet.found(curPoint))
1099  {
1100  pointSet.insert(curPoint);
1101  }
1102  }
1103  }
1104  }
1105 
1106  labelList corrPoints = pointSet.toc();
1107 
1108  labelListList corrPointFaces(corrPoints.size());
1109 
1110  forAll(corrPoints, pointI)
1111  {
1112  label curPoint = corrPoints[pointI];
1113 
1115 
1116  forAll(pFaces[curPoint], faceI)
1117  {
1118  label curFace = pFaces[curPoint][faceI];
1119 
1120  label index = eFaces.find(curFace);
1121 
1122  if (index != -1)
1123  {
1124  if (!faceSet.found(curFace))
1125  {
1126  faceSet.insert(curFace);
1127  }
1128  }
1129  }
1130 
1131  corrPointFaces[pointI] = faceSet.toc();
1132  }
1133 
1134  forAll(corrPoints, pointI)
1135  {
1136  label curPoint = corrPoints[pointI];
1137 
1138  scalar curDisp = 0;
1139 
1140  const labelList& curPointFaces = corrPointFaces[pointI];
1141 
1142  forAll(curPointFaces, faceI)
1143  {
1144  const face& curFace = faces[curPointFaces[faceI]];
1145 
1146  label ptInFace = curFace.which(curPoint);
1147  label next = curFace.nextLabel(ptInFace);
1148  label prev = curFace.prevLabel(ptInFace);
1149 
1150  vector a = points[next] - points[curPoint];
1151  vector b = points[prev] - points[curPoint];
1152  const vector& c = pointsDisplacementDir()[curPoint];
1153 
1154  curDisp += 2*sweptVolCorr[curPointFaces[faceI]]/((a^b)&c);
1155  }
1156 
1157  curDisp /= curPointFaces.size();
1158 
1159  displacement[curPoint] =
1160  curDisp*pointsDisplacementDir()[curPoint];
1161  }
1162  }
1163 
1164 
1165  for (const word& patchName : nonReflectingFreeSurfacePatches_)
1166  {
1167  label nonReflectingPatchID =
1168  aMesh().boundary().findPatchID(patchName);
1169 
1170  const labelList& pLabels =
1171  aMesh().boundary()[nonReflectingPatchID].pointLabels();
1172 
1173  const labelList& eFaces =
1174  aMesh().boundary()[nonReflectingPatchID].edgeFaces();
1175 
1176  labelList corrPoints = pLabels;
1177 
1178  labelListList corrPointFaces(corrPoints.size());
1179 
1180  forAll(corrPoints, pointI)
1181  {
1182  label curPoint = corrPoints[pointI];
1183 
1185 
1186  forAll(pFaces[curPoint], faceI)
1187  {
1188  label curFace = pFaces[curPoint][faceI];
1189 
1190  label index = eFaces.find(curFace);
1191 
1192  if (index != -1)
1193  {
1194  if (!faceSet.found(curFace))
1195  {
1196  faceSet.insert(curFace);
1197  }
1198  }
1199  }
1200 
1201  corrPointFaces[pointI] = faceSet.toc();
1202  }
1203 
1204 
1205  forAll(corrPoints, pointI)
1206  {
1207  label curPoint = corrPoints[pointI];
1208 
1209  scalar curDisp = 0;
1210 
1211  const labelList& curPointFaces = corrPointFaces[pointI];
1212 
1213  forAll(curPointFaces, faceI)
1214  {
1215  const face& curFace = faces[curPointFaces[faceI]];
1216 
1217  label ptInFace = curFace.which(curPoint);
1218  label next = curFace.nextLabel(ptInFace);
1219  label prev = curFace.prevLabel(ptInFace);
1220 
1221  label p0 = -1;
1222  label p1 = -1;
1223  label p2 = -1;
1224 
1225  if (corrPoints.find(next) == -1)
1226  {
1227  p0 = curPoint;
1228  p1 = next;
1229  p2 = curFace.nextLabel(curFace.which(next));
1230  }
1231  else
1232  {
1233  p0 = curFace.prevLabel(curFace.which(prev));
1234  p1 = prev;
1235  p2 = curPoint;
1236  }
1237 
1238  vector a0 = points[p1] - points[p0];
1239  vector b0 = points[p2] - points[p1];
1240  vector c0 = displacement[p1];
1241 
1242  scalar V0 = mag((a0^b0)&c0)/2;
1243 
1244  scalar DV = sweptVolCorr[curPointFaces[faceI]] - V0;
1245 
1246  if (corrPoints.find(prev) != -1)
1247  {
1248  vector a = points[curPoint] - points[prev];
1249  vector b =
1250  (points[next] + displacement[next])
1251  - points[curPoint];
1252  const vector& c = pointsDisplacementDir()[curPoint];
1253 
1254  curDisp += 2*DV/((a^b)&c);
1255  }
1256  else
1257  {
1258  vector a = points[curPoint]
1259  - (points[prev] + displacement[prev]);
1260  vector b = points[next] - points[curPoint];
1261  const vector& c = pointsDisplacementDir()[curPoint];
1262 
1263  curDisp += 2*DV/((a^b)&c);
1264  }
1265  }
1266 
1267  curDisp /= curPointFaces.size();
1268 
1269  displacement[curPoint] =
1270  curDisp*pointsDisplacementDir()[curPoint];
1271  }
1272  }
1273 }
1274 
1275 
1276 void Foam::interfaceTrackingFvMesh::correctContactLinePointNormals()
1277 {
1278  // Correct normals for contact line points
1279  // according to specified contact angle
1280 
1281  vectorField& N =
1282  const_cast<vectorField&>
1283  (
1285  );
1286 
1287  if (contactAnglePtr_ && correctContactLineNormals())
1288  {
1289  Info<< "Correcting contact line normals" << endl;
1290 
1291  vectorField oldPoints(aMesh().nPoints(), Zero);
1292 
1293  const labelList& meshPoints = aMesh().patch().meshPoints();
1294 
1295  forAll(oldPoints, ptI)
1296  {
1297  oldPoints[ptI] =
1298  mesh().oldPoints()[meshPoints[ptI]];
1299  }
1300 
1301 // # include "createTangentField.H"
1302  areaVectorField tangent
1303  (
1304  IOobject
1305  (
1306  "tangent",
1307  mesh().time().timeName(),
1308  mesh(),
1311  ),
1312  aMesh(),
1314  );
1315 
1316  if (Pstream::parRun())
1317  {
1318  const labelListList& edgeFaces = aMesh().patch().edgeFaces();
1319  const labelListList& pointEdges = aMesh().patch().pointEdges();
1320  const labelListList& pointFaces = aMesh().patch().pointFaces();
1321  const edgeList& edges = aMesh().edges();
1322 
1323  forAll(aMesh().boundary(), patchI)
1324  {
1325  if
1326  (
1327  aMesh().boundary()[patchI].type()
1328  == processorFaPatch::typeName
1329  )
1330  {
1331  const processorFaPatch& procPatch =
1332  refCast<const processorFaPatch>
1333  (
1334  aMesh().boundary()[patchI]
1335  );
1336 
1337  const labelList& patchPointLabels =
1338  procPatch.pointLabels();
1339 
1340  forAll(patchPointLabels, pointI)
1341  {
1342  label curPoint = patchPointLabels[pointI];
1343 
1344  // Check if processor point is boundary point
1345 
1346  label patchID = -1;
1347  label edgeID = -1;
1348 
1349  const labelList& curPointEdges = pointEdges[curPoint];
1350 
1351  forAll(curPointEdges, edgeI)
1352  {
1353  label curEdge = curPointEdges[edgeI];
1354 
1355  if (edgeFaces[curEdge].size() == 1)
1356  {
1357  forAll(aMesh().boundary(), pI)
1358  {
1359  const labelList& curEdges =
1360  aMesh().boundary()[pI];
1361 
1362  label index = curEdges.find(curEdge);
1363 
1364  if (index != -1)
1365  {
1366  if
1367  (
1368  aMesh().boundary()[pI].type()
1369  != processorFaPatch::typeName
1370  )
1371  {
1372  patchID = pI;
1373  edgeID = index;
1374  break;
1375  }
1376  }
1377  }
1378  }
1379  }
1380 
1381  if (patchID != -1)
1382  {
1383  label curEdge =
1384  aMesh().boundary()[patchID].start() + edgeID;
1385 
1386  vector t = edges[curEdge].vec(oldPoints);
1387  t /= mag(t) + SMALL;
1388 
1389  const labelList& curPointFaces =
1390  pointFaces[curPoint];
1391 
1392  forAll(curPointFaces, fI)
1393  {
1394  tangent.ref().field()[curPointFaces[fI]] = t;
1395  }
1396  }
1397  }
1398  }
1399  }
1400 
1401  tangent.correctBoundaryConditions();
1402  }
1403 
1404  forAll(aMesh().boundary(), patchI)
1405  {
1406  label ngbPolyPatchID =
1407  aMesh().boundary()[patchI].ngbPolyPatchIndex();
1408 
1409  if (ngbPolyPatchID != -1)
1410  {
1411  if
1412  (
1413  mesh().boundary()[ngbPolyPatchID].type()
1414  == wallFvPatch::typeName
1415  )
1416  {
1417  const scalar rotAngle = degToRad
1418  (
1419  gAverage
1420  (
1421  90
1422  - contactAnglePtr_->boundaryField()[patchI]
1423  )
1424  );
1425 
1426  vectorField ngbN
1427  (
1428  aMesh().boundary()[patchI].ngbPolyPatchPointNormals()
1429  );
1430 
1431  const labelList& patchPoints =
1432  aMesh().boundary()[patchI].pointLabels();
1433 
1434  vectorField pN(N, patchPoints);
1435 
1436  vectorField rotationAxis(ngbN^pN);
1437  rotationAxis /= mag(rotationAxis) + SMALL;
1438 
1439 
1440  // Calc rotation axis using edge vectors
1441 
1442  const edgeList& edges = aMesh().edges();
1443 
1444  const labelListList& pointEdges =
1445  aMesh().boundary()[patchI].pointEdges();
1446 
1447  forAll(pointEdges, pointI)
1448  {
1449  vector rotAx = Zero;
1450 
1451  forAll(pointEdges[pointI], eI)
1452  {
1453  label curEdge =
1454  aMesh().boundary()[patchI].start()
1455  + pointEdges[pointI][eI];
1456 
1457  vector e = edges[curEdge].vec(oldPoints);
1458 
1459  e *= (e&rotationAxis[pointI])
1460  /mag(e&rotationAxis[pointI]);
1461 
1462  e /= mag(e) + SMALL;
1463 
1464  rotAx += e;
1465  }
1466 
1467  if (pointEdges[pointI].size() == 1)
1468  {
1469  label curPoint = patchPoints[pointI];
1470 
1471  const labelListList& ptEdges =
1472  aMesh().patch().pointEdges();
1473  const labelList& curPointEdges =
1474  ptEdges[curPoint];
1475 
1476  label procPatchID = -1;
1477  label edgeID = -1;
1478 
1479  const labelListList& edgeFaces =
1480  aMesh().patch().edgeFaces();
1481 
1482  forAll(curPointEdges, edgeI)
1483  {
1484  label curEdge = curPointEdges[edgeI];
1485 
1486  if (edgeFaces[curEdge].size() == 1)
1487  {
1488  forAll(aMesh().boundary(), pI)
1489  {
1490  const labelList& curEdges =
1491  aMesh().boundary()[pI];
1492 
1493  label index =
1494  curEdges.find(curEdge);
1495 
1496  if (index != -1)
1497  {
1498  if
1499  (
1500  aMesh().boundary()[pI].type()
1501  == processorFaPatch::typeName
1502  )
1503  {
1504  procPatchID = pI;
1505  edgeID = index;
1506  break;
1507  }
1508  }
1509  }
1510  }
1511  }
1512 
1513  if (procPatchID != -1)
1514  {
1515  vector t =
1516  tangent.boundaryField()[procPatchID]
1517  .patchNeighbourField()()[edgeID];
1518 
1519  t *= (t&rotationAxis[pointI])
1520  /(mag(t&rotationAxis[pointI]) + SMALL);
1521 
1522  t /= mag(t) + SMALL;
1523 
1524  rotAx += t;
1525  }
1526  }
1527 
1528  rotationAxis[pointI] = rotAx/(mag(rotAx) + SMALL);
1529  }
1530 
1531  // Rodrigues' rotation formula
1532  ngbN = ngbN*cos(rotAngle)
1533  + rotationAxis*(rotationAxis & ngbN)*(1 - cos(rotAngle))
1534  + (rotationAxis^ngbN)*sin(rotAngle);
1535 
1536  // Info<< aMesh().boundary()[patchI].name() << endl;
1537  forAll(patchPoints, pointI)
1538  {
1539  N[patchPoints[pointI]] -=
1540  ngbN[pointI]*(ngbN[pointI]&N[patchPoints[pointI]]);
1541 
1542  N[patchPoints[pointI]] /=
1543  mag(N[patchPoints[pointI]]) + SMALL;
1544 
1545  // Info<< N[patchPoints[pointI]] << endl;
1546  }
1547  }
1548  }
1549  }
1550  }
1551 }
1552 
1553 
1554 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1555 
1556 Foam::interfaceTrackingFvMesh::interfaceTrackingFvMesh(const IOobject& io)
1557 :
1559  aMeshPtr_(new faMesh(*this)),
1560  fsPatchIndex_(-1),
1561  fixedFreeSurfacePatches_
1562  (
1563  motion().get<wordList>("fixedFreeSurfacePatches")
1564  ),
1565  nonReflectingFreeSurfacePatches_(),
1566  pointNormalsCorrectionPatches_
1567  (
1568  motion().get<wordList>("pointNormalsCorrectionPatches")
1569  ),
1570  normalMotionDir_
1571  (
1572  motion().get<bool>("normalMotionDir")
1573  ),
1574  motionDir_(Zero),
1575  smoothing_(motion().getOrDefault("smoothing", false)),
1576  pureFreeSurface_(motion().getOrDefault("pureFreeSurface", true)),
1577  rigidFreeSurface_(motion().getOrDefault("rigidFreeSurface", false)),
1578  correctContactLineNormals_
1579  (
1580  motion().getOrDefault("correctContactLineNormals", false)
1581  ),
1582  sigma0_("zero", dimForce/dimLength/dimDensity, Zero),
1583  rho_("one", dimDensity, 1.0),
1584  timeIndex_(-1),
1585  UsPtr_(nullptr),
1586  controlPointsPtr_(nullptr),
1587  motionPointsMaskPtr_(nullptr),
1588  pointsDisplacementDirPtr_(nullptr),
1589  facesDisplacementDirPtr_(nullptr),
1590  fsNetPhiPtr_(nullptr),
1591  phisPtr_(nullptr),
1592  surfactConcPtr_(nullptr),
1593  bulkSurfactConcPtr_(nullptr),
1594  surfaceTensionPtr_(nullptr),
1595  surfactantPtr_(nullptr),
1596  contactAnglePtr_(nullptr)
1597 {
1598  initializeData();
1599 }
1600 
1601 
1602 Foam::interfaceTrackingFvMesh::interfaceTrackingFvMesh
1604  const IOobject& io,
1605  pointField&& points,
1606  faceList&& faces,
1607  labelList&& allOwner,
1608  labelList&& allNeighbour,
1609  const bool syncPar
1610 )
1611 :
1613  (
1614  io,
1615  std::move(points),
1616  std::move(faces),
1617  std::move(allOwner),
1618  std::move(allNeighbour),
1619  syncPar
1620  ),
1621  aMeshPtr_(new faMesh(*this)),
1622  fsPatchIndex_(-1),
1623  fixedFreeSurfacePatches_
1624  (
1625  motion().get<wordList>("fixedFreeSurfacePatches")
1626  ),
1627  nonReflectingFreeSurfacePatches_(),
1628  pointNormalsCorrectionPatches_
1629  (
1630  motion().get<wordList>("pointNormalsCorrectionPatches")
1631  ),
1632  normalMotionDir_
1633  (
1634  motion().get<bool>("normalMotionDir")
1635  ),
1636  motionDir_(Zero),
1637  smoothing_(motion().getOrDefault("smoothing", false)),
1638  pureFreeSurface_(motion().getOrDefault("pureFreeSurface", true)),
1639  sigma0_("zero", dimForce/dimLength/dimDensity, Zero),
1640  rho_("one", dimDensity, 1.0),
1641  timeIndex_(-1),
1642  UsPtr_(nullptr),
1643  controlPointsPtr_(nullptr),
1644  motionPointsMaskPtr_(nullptr),
1645  pointsDisplacementDirPtr_(nullptr),
1646  facesDisplacementDirPtr_(nullptr),
1647  fsNetPhiPtr_(nullptr),
1648  phisPtr_(nullptr),
1649  surfactConcPtr_(nullptr),
1650  bulkSurfactConcPtr_(nullptr),
1651  surfaceTensionPtr_(nullptr),
1652  surfactantPtr_(nullptr),
1653  contactAnglePtr_(nullptr)
1654 {
1655  initializeData();
1656 }
1657 
1658 
1659 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1660 
1662 {
1663  deleteDemandDrivenData(UsPtr_);
1664  deleteDemandDrivenData(controlPointsPtr_);
1665  deleteDemandDrivenData(motionPointsMaskPtr_);
1666  deleteDemandDrivenData(pointsDisplacementDirPtr_);
1667  deleteDemandDrivenData(facesDisplacementDirPtr_);
1668  deleteDemandDrivenData(fsNetPhiPtr_);
1669  deleteDemandDrivenData(phisPtr_);
1670  deleteDemandDrivenData(surfactConcPtr_);
1671  deleteDemandDrivenData(bulkSurfactConcPtr_);
1672  deleteDemandDrivenData(surfaceTensionPtr_);
1673  deleteDemandDrivenData(surfactantPtr_);
1674  deleteDemandDrivenData(contactAnglePtr_);
1675 }
1676 
1677 
1678 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1679 
1681 {
1682  if (!UsPtr_)
1683  {
1684  makeUs();
1685  }
1686 
1687  return *UsPtr_;
1688 }
1689 
1690 
1692 {
1693  if (!UsPtr_)
1694  {
1695  makeUs();
1696  }
1697 
1698  return *UsPtr_;
1699 }
1700 
1701 
1703 {
1704  if (!fsNetPhiPtr_)
1705  {
1706  makeFsNetPhi();
1707  }
1708 
1709  return *fsNetPhiPtr_;
1710 }
1711 
1712 
1714 {
1715  if (!fsNetPhiPtr_)
1716  {
1717  makeFsNetPhi();
1718  }
1719 
1720  return *fsNetPhiPtr_;
1721 }
1722 
1723 
1725 {
1726  forAll(Us().boundaryField(), patchI)
1727  {
1728  if
1729  (
1730  Us().boundaryField()[patchI].type()
1731  == calculatedFaPatchVectorField::typeName
1732  )
1733  {
1734  vectorField& pUs = Us().boundaryFieldRef()[patchI];
1735 
1736  pUs = Us().boundaryField()[patchI].patchInternalField();
1737 
1738  label ngbPolyPatchID =
1739  aMesh().boundary()[patchI].ngbPolyPatchIndex();
1740 
1741  if (ngbPolyPatchID != -1)
1742  {
1743  if
1744  (
1745  (
1746  U().boundaryField()[ngbPolyPatchID].type()
1747  == slipFvPatchVectorField::typeName
1748  )
1749  ||
1750  (
1751  U().boundaryField()[ngbPolyPatchID].type()
1752  == symmetryFvPatchVectorField::typeName
1753  )
1754  )
1755  {
1756  vectorField N
1757  (
1758  aMesh().boundary()[patchI].ngbPolyPatchFaceNormals()
1759  );
1760 
1761  pUs -= N*(N&pUs);
1762  }
1763  }
1764  }
1765  }
1766 
1767  Us().correctBoundaryConditions();
1768 }
1769 
1770 
1772 {
1773  // Info<< "Update Us" << endl;
1774 
1775  Us().ref().field() = U().boundaryField()[fsPatchIndex()];
1776 
1777  // // Correct normal component of free-surface velocity
1778  // const vectorField& nA = aMesh().faceAreaNormals().internalField();
1779  // vectorField UnFs = nA*phi().boundaryField()[fsPatchIndex()]
1780  // /mesh().boundary()[fsPatchIndex()].magSf();
1781  // Us().ref().field() += UnFs - nA*(nA&Us().internalField());
1782 
1783  correctUsBoundaryConditions();
1784 }
1785 
1786 
1788 {
1789  return *getObjectPtr<const volVectorField>("U");
1790 }
1791 
1792 
1794 {
1795  return *getObjectPtr<const volScalarField>("p");
1796 }
1797 
1798 
1800 {
1801  return *getObjectPtr<const surfaceScalarField>("phi");
1802 }
1803 
1804 
1807 {
1808  auto tSnGradU = tmp<vectorField>::New(aMesh().nFaces(), Zero);
1809  auto& SnGradU = tSnGradU.ref();
1810 
1811  const vectorField& nA = aMesh().faceAreaNormals().internalField();
1812 
1813  areaScalarField divUs
1814  (
1815  fac::div(Us())
1816  - aMesh().faceCurvatures()*(aMesh().faceAreaNormals()&Us())
1817  );
1818 
1819  areaTensorField gradUs(fac::grad(Us()));
1820 
1821  // Remove component of gradient normal to surface (area)
1822  const areaVectorField& n = aMesh().faceAreaNormals();
1823  gradUs -= n*(n & gradUs);
1824  gradUs.correctBoundaryConditions();
1825 
1826  const turbulenceModel& turbulence =
1827  mesh().lookupObject<turbulenceModel>("turbulenceProperties");
1828 
1829  scalarField nu(turbulence.nuEff(fsPatchIndex()));
1830 
1831  vectorField tangentialSurfaceTensionForce(nA.size(), Zero);
1832 
1833  if (!pureFreeSurface() && max(nu) > SMALL)
1834  {
1835  tangentialSurfaceTensionForce =
1836  surfaceTensionGrad()().internalField();
1837  }
1838 
1839  SnGradU =
1840  tangentialSurfaceTensionForce/(nu + SMALL)
1841  - nA*divUs.internalField()
1842  - (gradUs.internalField()&nA);
1843 
1844  return tSnGradU;
1845 }
1846 
1847 
1850 {
1851  auto tSnGradUn = tmp<scalarField>::New(aMesh().nFaces(), Zero);
1852  auto& SnGradUn = tSnGradUn.ref();
1853 
1854  areaScalarField divUs
1855  (
1856  fac::div(Us())
1857  - aMesh().faceCurvatures()*(aMesh().faceAreaNormals()&Us())
1858  );
1859 
1860  SnGradUn = -divUs.internalField();
1861 
1862  return tSnGradUn;
1863 }
1864 
1865 
1868 {
1869  auto tPressureJump = tmp<scalarField>::New(aMesh().nFaces(), Zero);
1870  auto& pressureJump = tPressureJump.ref();
1871 
1873 
1875  meshObjects::gravity::New(mesh().time());
1876 
1877  const turbulenceModel& turbulence =
1878  mesh().lookupObject<turbulenceModel>("turbulenceProperties");
1879 
1880  scalarField nu(turbulence.nuEff(fsPatchIndex()));
1881 
1882  pressureJump =
1883  - (g.value() & mesh().Cf().boundaryField()[fsPatchIndex()])
1884  + 2.0*nu*freeSurfaceSnGradUn();
1885 
1886  if (pureFreeSurface())
1887  {
1888  pressureJump -= sigma().value()*K;
1889  }
1890  else
1891  {
1892  pressureJump -= surfaceTension().internalField()*K;
1893  }
1894 
1895  return tPressureJump;
1896 }
1897 
1898 
1900 {
1901  if (!controlPointsPtr_)
1902  {
1903  makeControlPoints();
1904  }
1905 
1906  return *controlPointsPtr_;
1907 }
1908 
1909 
1911 {
1912  if (!motionPointsMaskPtr_)
1913  {
1914  makeMotionPointsMask();
1915  }
1916 
1917  return *motionPointsMaskPtr_;
1918 }
1919 
1920 
1922 {
1923  if (!pointsDisplacementDirPtr_)
1924  {
1925  makeDirections();
1926  }
1927 
1928  return *pointsDisplacementDirPtr_;
1929 }
1930 
1931 
1933 {
1934  if (!facesDisplacementDirPtr_)
1935  {
1936  makeDirections();
1937  }
1938 
1939  return *facesDisplacementDirPtr_;
1940 }
1941 
1942 
1944 {
1945  if (!phisPtr_)
1946  {
1947  makePhis();
1948  }
1949 
1950  return *phisPtr_;
1951 }
1952 
1955 {
1956  if (!surfactConcPtr_)
1957  {
1958  makeSurfactConc();
1959  }
1960 
1961  return *surfactConcPtr_;
1962 }
1963 
1964 
1965 const Foam::areaScalarField&
1967 {
1968  if (!surfactConcPtr_)
1969  {
1970  makeSurfactConc();
1971  }
1972 
1973  return *surfactConcPtr_;
1974 }
1975 
1976 
1979 {
1980  if (!bulkSurfactConcPtr_)
1981  {
1982  makeBulkSurfactConc();
1983  }
1984 
1985  return *bulkSurfactConcPtr_;
1986 }
1987 
1988 
1989 const Foam::volScalarField&
1991 {
1992  if (!bulkSurfactConcPtr_)
1993  {
1994  makeBulkSurfactConc();
1995  }
1996 
1997  return *bulkSurfactConcPtr_;
1998 }
1999 
2000 
2003 {
2004  if (!surfaceTensionPtr_)
2005  {
2006  makeSurfaceTension();
2007  }
2008 
2009  return *surfaceTensionPtr_;
2010 }
2011 
2012 
2013 const Foam::areaScalarField&
2015 {
2016  if (!surfaceTensionPtr_)
2017  {
2018  makeSurfaceTension();
2019  }
2020 
2021  return *surfaceTensionPtr_;
2022 }
2023 
2024 
2027 {
2028  if (!surfactantPtr_)
2029  {
2030  makeSurfactant();
2031  }
2032 
2033  return *surfactantPtr_;
2034 }
2035 
2036 
2039 {
2040  auto tgrad = tmp<areaVectorField>::New
2041  (
2042  IOobject
2043  (
2044  "surfaceTensionGrad",
2045  mesh().time().timeName(),
2046  mesh(),
2049  ),
2050  fac::grad(surfaceTension())
2051  // (-fac::grad(surfactantConcentration()/rho_)*
2052  // surfactant().surfactR()*surfactant().surfactT()/
2053  // (1.0 - surfactantConcentration()/
2054  // surfactant().surfactSaturatedConc()))()
2055  );
2056  auto& grad = tgrad.ref();
2057 
2058  // Remove component of gradient normal to surface (area)
2059  const areaVectorField& n = aMesh().faceAreaNormals();
2060  grad -= n*(n & grad);
2061  grad.correctBoundaryConditions();
2062 
2063  return tgrad;
2064 }
2065 
2066 
2068 {
2069  if (timeIndex_ != mesh().time().timeIndex())
2070  {
2071  if (smoothing_ && !rigidFreeSurface_)
2072  {
2073  smoothFreeSurfaceMesh();
2074  clearControlPoints();
2075  }
2076 
2077  updateDisplacementDirections();
2078 
2079  updateProperties();
2080 
2081  Info<< "Maximal capillary Courant number: "
2082  << maxCourantNumber() << endl;
2083 
2085 
2086  Info<< "Free surface curvature: min = " << gMin(K)
2087  << ", max = " << gMax(K) << ", average = " << gAverage(K) << nl;
2088 
2089  timeIndex_ = mesh().time().timeIndex();
2090  }
2091 
2092  if (!rigidFreeSurface_)
2093  {
2094  // This is currently relaltive flux
2095  scalarField sweptVolCorr =
2096  phi().boundaryField()[fsPatchIndex()];
2097 
2098  // Info<< "Free surface flux: sum local = "
2099  // << gSum(mag(sweptVolCorr))
2100  // << ", global = " << gSum(sweptVolCorr) << endl;
2101 
2102  // if (mesh().moving())
2103  // {
2104  // sweptVolCorr -=
2105  // fvc::meshPhi(U())().boundaryField()[fsPatchIndex()];
2106  // }
2107 
2108  Info<< "Free surface continuity error : sum local = "
2109  << gSum(mag(sweptVolCorr)) << ", global = " << gSum(sweptVolCorr)
2110  << endl;
2111 
2112  // For postprocessing
2113  fsNetPhi().ref().field() = sweptVolCorr;
2114 
2115  word ddtScheme
2116  (
2117  mesh().ddtScheme
2118  (
2119  "ddt(" + U().name() + ')'
2120  )
2121  );
2122 
2123  if
2124  (
2125  ddtScheme
2127  )
2128  {
2129  sweptVolCorr *= (1.0/2.0)*mesh().time().deltaT().value();
2130  }
2131  else if
2132  (
2133  ddtScheme
2135  )
2136  {
2137  sweptVolCorr *= mesh().time().deltaT().value();
2138  }
2139  else if
2140  (
2141  ddtScheme
2143  )
2144  {
2145  if (mesh().time().timeIndex() == 1)
2146  {
2147  sweptVolCorr *= mesh().time().deltaT().value();
2148  }
2149  else
2150  {
2151  sweptVolCorr *= (2.0/3.0)*mesh().time().deltaT().value();
2152  }
2153  }
2154  else
2155  {
2157  << "Unsupported temporal differencing scheme : "
2158  << ddtScheme << nl
2159  << abort(FatalError);
2160  }
2161 
2162  const scalarField& Sf = aMesh().S();
2163  const vectorField& Nf = aMesh().faceAreaNormals().internalField();
2164 
2165  scalarField deltaHf
2166  (
2167  sweptVolCorr/(Sf*(Nf & facesDisplacementDir()))
2168  );
2169 
2170  for (const word& patchName : fixedFreeSurfacePatches_)
2171  {
2172  label fixedPatchID =
2173  aMesh().boundary().findPatchID(patchName);
2174 
2175  if (fixedPatchID == -1)
2176  {
2178  << "Wrong faPatch name '" << patchName
2179  << "'in the fixedFreeSurfacePatches list"
2180  << " defined in dynamicMeshDict dictionary"
2181  << abort(FatalError);
2182  }
2183 
2184  const labelList& eFaces =
2185  aMesh().boundary()[fixedPatchID].edgeFaces();
2186 
2187  forAll(eFaces, edgeI)
2188  {
2189  deltaHf[eFaces[edgeI]] *= 2.0;
2190  }
2191  }
2192 
2193  controlPoints() += facesDisplacementDir()*deltaHf;
2194 
2195  pointField displacement(pointDisplacement());
2196  correctPointDisplacement(sweptVolCorr, displacement);
2197 
2198  velocityMotionSolver& vMotion =
2199  refCast<velocityMotionSolver>
2200  (
2201  const_cast<motionSolver&>(motion())
2202  );
2203 
2204  pointVectorField& pointMotionU = vMotion.pointMotionU();
2205  pointMotionU.primitiveFieldRef() = Zero;
2206 
2207  fixedValuePointPatchVectorField& fsPatchPointMeshU =
2208  refCast<fixedValuePointPatchVectorField>
2209  (
2210  const_cast<pointPatchVectorField&>
2211  (
2212  pointMotionU.boundaryField()[fsPatchIndex()]
2213  )
2214  );
2215 
2216  fsPatchPointMeshU ==
2217  displacement/mesh().time().deltaT().value();
2218 
2220 
2221  correctContactLinePointNormals();
2222  }
2223  else
2224  {
2225  vectorField displacement
2226  (
2227  mesh().boundaryMesh()[fsPatchIndex()].nPoints(),
2228  Zero
2229  );
2230 
2231  velocityMotionSolver& vMotion =
2232  refCast<velocityMotionSolver>
2233  (
2234  const_cast<motionSolver&>(motion())
2235  );
2236 
2237  pointVectorField& pointMotionU = vMotion.pointMotionU();
2238  pointMotionU.primitiveFieldRef() = Zero;
2239 
2240  fixedValuePointPatchVectorField& fsPatchPointMeshU =
2241  refCast<fixedValuePointPatchVectorField>
2242  (
2243  const_cast<pointPatchVectorField&>
2244  (
2245  pointMotionU.boundaryField()[fsPatchIndex()]
2246  )
2247  );
2248 
2249  fsPatchPointMeshU ==
2250  displacement/mesh().time().deltaT().value();
2251 
2253  }
2254 
2255  updateUs();
2256 
2257  updateSurfactantConcentration();
2258 
2259  return true;
2260 }
2261 
2262 
2264 {
2265  // Write patch and points into VTK
2266  OFstream mps(mesh().time().timePath()/"freeSurface.vtk");
2267 
2268  const vectorField& points = aMesh().patch().points();
2269  const IndirectList<face>& faces = aMesh().patch();
2270 
2271  mps << "# vtk DataFile Version 2.0" << nl
2272  << mesh().time().timePath()/"freeSurface.vtk" << nl
2273  << "ASCII" << nl
2274  << "DATASET POLYDATA" << nl
2275  << "POINTS " << points.size() << " float" << nl;
2276 
2277  // Write points
2278  List<float> mlpBuffer(3*points.size());
2279 
2280  label counter = 0;
2281  forAll(points, i)
2282  {
2283  mlpBuffer[counter++] = float(points[i].x());
2284  mlpBuffer[counter++] = float(points[i].y());
2285  mlpBuffer[counter++] = float(points[i].z());
2286  }
2287 
2288  forAll(mlpBuffer, i)
2289  {
2290  mps << mlpBuffer[i] << ' ';
2291 
2292  if (i > 0 && (i % 10) == 0)
2293  {
2294  mps << nl;
2295  }
2296  }
2297 
2298  // Write faces
2299  label nFaceVerts = 0;
2300 
2301  forAll(faces, faceI)
2302  {
2303  nFaceVerts += faces[faceI].size() + 1;
2304  }
2305  labelList mlfBuffer(nFaceVerts);
2306 
2307  counter = 0;
2308  forAll(faces, faceI)
2309  {
2310  const face& f = faces[faceI];
2311 
2312  mlfBuffer[counter++] = f.size();
2313 
2314  forAll(f, fpI)
2315  {
2316  mlfBuffer[counter++] = f[fpI];
2317  }
2318  }
2319  mps << nl;
2320 
2321  mps << "POLYGONS " << faces.size() << ' ' << nFaceVerts << endl;
2322 
2323  forAll(mlfBuffer, i)
2324  {
2325  mps << mlfBuffer[i] << ' ';
2326 
2327  if (i > 0 && (i % 10) == 0)
2328  {
2329  mps << nl;
2330  }
2331  }
2332  mps << nl;
2333 
2334  // aMesh().patch().writeVTK
2335  // (
2336  // mesh().time().timePath()/"freeSurface",
2337  // aMesh().patch(),
2338  // aMesh().patch().points()
2339  // );
2340 }
2341 
2342 
2344 {
2345  // Write control points into VTK
2346  fileName name(mesh().time().timePath()/"freeSurfaceControlPoints.vtk");
2347  OFstream mps(name);
2348 
2349  Info<< "Writing free surface control point to " << name << endl;
2350 
2351  mps << "# vtk DataFile Version 2.0" << nl
2352  << name << nl
2353  << "ASCII" << nl
2354  << "DATASET POLYDATA" << nl
2355  << "POINTS " << controlPoints().size() << " float" << nl;
2356 
2357  forAll(controlPoints(), pointI)
2358  {
2359  mps << controlPoints()[pointI].x() << ' '
2360  << controlPoints()[pointI].y() << ' '
2361  << controlPoints()[pointI].z() << nl;
2362  }
2363 
2364  // Write vertices
2365  mps << "VERTICES " << controlPoints().size() << ' '
2366  << controlPoints().size()*2 << nl;
2367 
2368  forAll(controlPoints(), pointI)
2369  {
2370  mps << 1 << ' ' << pointI << nl;
2371  }
2372 }
2373 
2374 
2375 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
volFields.H
Foam::vectorIOField
IOField< vector > vectorIOField
vectorField with IO.
Definition: vectorIOField.H:44
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
backwardDdtScheme.H
Foam::faMesh::Le
const edgeVectorField & Le() const
Return edge length vectors.
Definition: faMesh.C:1036
Foam::interfaceTrackingFvMesh::aMesh
faMesh & aMesh()
Return reference to finite area mesh.
Definition: interfaceTrackingFvMesh.H:284
Foam::PrimitivePatch::pointFaces
const labelListList & pointFaces() const
Return point-face addressing.
Definition: PrimitivePatch.C:281
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::interfaceTrackingFvMesh::bulkSurfactantConcentration
volScalarField & bulkSurfactantConcentration()
Return volume surfactant concentration field.
Definition: interfaceTrackingFvMesh.C:1978
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatch.C:242
wedgeFaPatchFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::faMesh::faces
const faceList & faces() const
Return faces.
Definition: faMesh.C:965
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::TimeState::deltaT
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:54
Foam::HashTable::toc
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:121
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::interfaceTrackingFvMesh::surfactantConcentration
areaScalarField & surfactantConcentration()
Return free-surface surfactant concentration field.
Definition: interfaceTrackingFvMesh.C:1954
motionSolver.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::faMesh::patch
const indirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMesh.C:915
slipFvPatchFields.H
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::fam::laplacian
tmp< faMatrix< Type > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famLaplacian.C:49
Foam::dimDensity
const dimensionSet dimDensity
Foam::dynamicMotionSolverFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: dynamicMotionSolverFvMesh.C:93
fvcMeshPhi.H
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::faMesh::correctPatchPointNormals
bool correctPatchPointNormals(const label patchID) const
Whether point normals should be corrected for a patch.
Definition: faMesh.C:1276
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:415
turbulentTransportModel.H
Foam::interfaceTrackingFvMesh
The interfaceTrackingFvMesh.
Definition: interfaceTrackingFvMesh.H:57
gravityMeshObject.H
Foam::dynamicFvMesh
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:78
wallFvPatch.H
Foam::interfaceTrackingFvMesh::writeVTK
void writeVTK() const
Write VTK freeSurface mesh.
Definition: interfaceTrackingFvMesh.C:2263
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
unitConversion.H
Unit conversion functions.
Foam::faMesh::edges
const edgeList & edges() const
Return edges.
Definition: faMesh.C:959
Foam::faceSet
A list of face labels.
Definition: faceSet.H:51
Foam::dimMoles
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:56
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::faMesh::pointAreaNormals
const vectorField & pointAreaNormals() const
Return point area normals.
Definition: faMesh.C:1153
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
faScalarMatrix
Template specialisation for scalar faMatrix.
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:57
Foam::GeometricField::internalField
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Definition: GeometricFieldI.H:43
Foam::dimForce
const dimensionSet dimForce
Foam::velocityMotionSolver
Virtual base class for velocity motion solver.
Definition: velocityMotionSolver.H:58
Foam::DynamicID< polyBoundaryMesh >
Foam::dynamicMotionSolverFvMesh::motion
const motionSolver & motion() const
Return the motionSolver.
Definition: dynamicMotionSolverFvMesh.C:87
Foam::HashSet< label, Hash< label > >
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
Foam::pointPatchField< vector >
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
polyPatchID.H
Foam::fam::ddt
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:48
edgeID
label edgeID
Definition: boundaryProcessorFaPatchPoints.H:6
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::interfaceTrackingFvMesh::freeSurfacePressureJump
tmp< scalarField > freeSurfacePressureJump()
Return free surface pressure jump.
Definition: interfaceTrackingFvMesh.C:1867
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
aMesh
faMesh aMesh(mesh)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::areaVectorField
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:56
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::flush
Ostream & flush(Ostream &os)
Flush stream.
Definition: Ostream.H:342
Foam::processorFaPatch
Processor patch.
Definition: processorFaPatch.H:55
Foam::fam::Sp
tmp< faMatrix< Type > > Sp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:86
slipFaPatchFields.H
Foam::interfaceTrackingFvMesh::updateUs
void updateUs()
Update free-surface velocity field.
Definition: interfaceTrackingFvMesh.C:1771
Foam::topoSet::found
virtual bool found(const label id) const
Has the given index?
Definition: topoSet.C:511
Foam::UniformDimensionedField< vector >
velocityLaplacianFvMotionSolver.H
Foam::interfaceTrackingFvMesh::surfaceTensionGrad
tmp< areaVectorField > surfaceTensionGrad()
Return surface tension gradient.
Definition: interfaceTrackingFvMesh.C:2038
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::face::which
label which(const label pointLabel) const
Find local index on face for the point label,.
Definition: faceI.H:138
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::interfaceTrackingFvMesh::p
const volScalarField & p() const
Return constant reference to pressure field.
Definition: interfaceTrackingFvMesh.C:1793
Foam::faMesh::S
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:1081
Foam::fv::CrankNicolsonDdtScheme
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Definition: CrankNicolsonDdtScheme.H:115
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:360
EulerDdtScheme.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::Time::timePath
fileName timePath() const
Return current time path.
Definition: Time.H:375
Foam::faMesh::faceAreaNormals
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:1131
Foam::edgeInterpolation::lPN
const edgeScalarField & lPN() const
Return reference to PN geodesic distance.
Definition: edgeInterpolation.C:82
Foam::interfaceTrackingFvMesh::pointsDisplacementDir
vectorField & pointsDisplacementDir()
Return reference to point displacement direction field.
Definition: interfaceTrackingFvMesh.C:1921
fixedValueFaPatchFields.H
Foam::PrimitivePatch::pointEdges
const labelListList & pointEdges() const
Return point-edge addressing.
Definition: PrimitivePatch.C:268
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
field
rDeltaTY field()
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::edgeScalarField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:52
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::polyMesh::oldPoints
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1088
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fac::edgeIntegrate
tmp< GeometricField< Type, faPatchField, areaMesh > > edgeIntegrate(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facEdgeIntegrate.C:47
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::faBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name, return -1 if not found.
Definition: faBoundaryMesh.C:290
CrankNicolsonDdtScheme.H
interfaceTrackingFvMesh.H
Foam::PrimitivePatch::points
const Field< point_type > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:305
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::interfaceTrackingFvMesh::U
const volVectorField & U() const
Return constant reference to velocity field.
Definition: interfaceTrackingFvMesh.C:1787
Foam::dimensioned< vector >
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:297
Foam::interfaceTrackingFvMesh::fsNetPhi
const areaScalarField & fsNetPhi() const
Return free-surface net flux.
Definition: interfaceTrackingFvMesh.C:1713
Foam::interfaceTrackingFvMesh::freeSurfaceSnGradU
tmp< vectorField > freeSurfaceSnGradU()
Return free surface normal derivative of velocity.
Definition: interfaceTrackingFvMesh.C:1806
Foam::faMesh::points
const pointField & points() const
Return mesh points.
Definition: faMesh.C:953
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::surfactantProperties
Definition: surfactantProperties.H:50
Foam::areaScalarField
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:53
solveBulkSurfactant.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::interfaceTrackingFvMesh::freeSurfaceSnGradUn
tmp< scalarField > freeSurfaceSnGradUn()
Return free surface normal derivative of normal velocity comp.
Definition: interfaceTrackingFvMesh.C:1849
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
CoNum
scalar CoNum
Definition: compressibleCourantNo.H:31
Foam::interfaceTrackingFvMesh::surfaceTension
areaScalarField & surfaceTension()
Return surface tension field.
Definition: interfaceTrackingFvMesh.C:2002
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:339
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:87
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:483
Foam::linearEdgeInterpolate
tmp< GeometricField< Type, faePatchField, edgeMesh > > linearEdgeInterpolate(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: linearEdgeInterpolation.H:108
Foam::faPatch::pointLabels
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:138
Foam::interfaceTrackingFvMesh::motionPointsMask
labelList & motionPointsMask()
Return reference to motion points mask field.
Definition: interfaceTrackingFvMesh.C:1910
U
U
Definition: pEqn.H:72
Foam::face::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:150
Foam::constant::atomic::a0
const dimensionedScalar a0
Bohr radius: default SI units: [m].
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::motionSolver
Virtual base class for mesh motion solver.
Definition: motionSolver.H:58
CGAL::indexedCell
An indexed form of CGAL::Triangulation_cell_base_3<K> used to keep track of the Delaunay cells (tets)...
Definition: indexedCell.H:58
Foam::faMesh::faceCurvatures
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition: faMesh.C:1169
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::pointSet
A set of point labels.
Definition: pointSet.H:51
Foam::fv::backwardDdtScheme
Second-order backward-differencing ddt using the current and two previous time-step values.
Definition: backwardDdtScheme.H:61
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List< bool >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::interfaceTrackingFvMesh::writeVTKControlPoints
void writeVTKControlPoints()
Write VTK freeSurface control points.
Definition: interfaceTrackingFvMesh.C:2343
Us
Us
Definition: createFaFields.H:51
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
symmetryFvPatchFields.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::interfaceTrackingFvMesh::facesDisplacementDir
vectorField & facesDisplacementDir()
Return reference to control points displacement direction field.
Definition: interfaceTrackingFvMesh.C:1932
Foam::faMesh::boundary
const faBoundaryMesh & boundary() const
Return constant reference to boundary mesh.
Definition: faMesh.C:1019
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
Foam::interfaceTrackingFvMesh::correctUsBoundaryConditions
void correctUsBoundaryConditions()
Correct surface velocity boundary conditions.
Definition: interfaceTrackingFvMesh.C:1724
bool
bool
Definition: EEqn.H:20
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::velocityMotionSolver::pointMotionU
pointVectorField & pointMotionU()
Return reference to the point motion velocity field.
Definition: velocityMotionSolver.H:102
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::interfaceTrackingFvMesh::phi
const surfaceScalarField & phi() const
Return constant reference to velocity field.
Definition: interfaceTrackingFvMesh.C:1799
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:77
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
twoDPointCorrector.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Cb
CGAL::indexedCell< K > Cb
Definition: CGALTriangulation3Ddefs.H:49
startTime
Foam::label startTime
Definition: checkTimeOptions.H:1
Foam::interfaceTrackingFvMesh::controlPoints
vectorField & controlPoints()
Return control points.
Definition: interfaceTrackingFvMesh.C:1899
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::interfaceTrackingFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: interfaceTrackingFvMesh.C:2067
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
N
const Vector< label > N(dict.get< Vector< label >>("N"))
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:36
Foam::interfaceTrackingFvMesh::surfactant
const surfactantProperties & surfactant() const
Return surfactant properties.
Definition: interfaceTrackingFvMesh.C:2026
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:310
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:352
wedgeFaPatch.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fv::EulerDdtScheme
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Definition: EulerDdtScheme.H:60
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::interfaceTrackingFvMesh::Us
areaVectorField & Us()
Return free-surface velocity field.
Definition: interfaceTrackingFvMesh.C:1680
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::interfaceTrackingFvMesh::Phis
edgeScalarField & Phis()
Return free-surface fluid flux field.
Definition: interfaceTrackingFvMesh.C:1943
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::interfaceTrackingFvMesh::~interfaceTrackingFvMesh
~interfaceTrackingFvMesh()
Destructor.
Definition: interfaceTrackingFvMesh.C:1661
Foam::dynamicMotionSolverFvMesh
The dynamicMotionSolverFvMesh.
Definition: dynamicMotionSolverFvMesh.H:52
Foam::fam::div
tmp< faMatrix< Type > > div(const edgeScalarField &flux, const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: famDiv.C:48
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::faMesh::areaCentres
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition: faMesh.C:1058
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265