chombo-discharge
Loading...
Searching...
No Matches
CD_AmrMeshImplem.H
Go to the documentation of this file.
1/* chombo-discharge
2 * Copyright © 2021 SINTEF Energy Research.
3 * Please refer to Copyright.txt and LICENSE in the chombo-discharge root directory.
4 */
5
12#ifndef CD_AmrMeshImplem_H
13#define CD_AmrMeshImplem_H
14
15// Our includes
16#include <CD_AmrMesh.H>
17#include <CD_ParticleOps.H>
18#include <CD_NamespaceHeader.H>
19
20template <typename T>
21void
23 const EBAMRData<T>& a_src,
25 const CopyStrategy& a_fromRegion) const noexcept
26{
27 CH_TIME("AmrMesh::copyData(EBAMRData, simple)");
28 if (m_verbosity > 5) {
29 pout() << "AmrMesh::copyData(EBAMRData, simple)" << endl;
30 }
31
32 const int nComp = a_dst[0]->nComp();
33 const Interval dstComps = Interval(0, nComp - 1);
34 const Interval srcComps = Interval(0, nComp - 1);
35
36 const std::string toRealm = a_dst.getRealm();
37 const std::string fromRealm = a_src.getRealm();
38
39 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
42
44 }
45}
46
47template <typename T>
48void
50 const EBAMRData<T>& a_src,
51 const Interval& a_dstComps,
52 const Interval& a_srcComps,
54 const CopyStrategy& a_fromRegion) const noexcept
55{
56 CH_TIME("AmrMesh::copyData(EBAMRData, full)");
57 if (m_verbosity > 5) {
58 pout() << "AmrMesh::copyData(EBAMRData, full)" << endl;
59 }
60
61 const std::string toRealm = a_dst.getRealm();
62 const std::string fromRealm = a_src.getRealm();
63
64 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
67
69 }
70}
71
72template <typename T>
73void
75 const LevelData<T>& a_src,
76 const int a_level,
80 const CopyStrategy& a_fromRegion) const noexcept
81{
82 CH_TIME("AmrMesh::copyData(LD<T>_full)");
83 if (m_verbosity > 5) {
84 pout() << "AmrMesh::copyData(LD<T>_full)" << endl;
85 }
86
87 const Interval dstComps = Interval(0, a_dst.nComp() - 1);
88 const Interval srcComps = Interval(0, a_src.nComp() - 1);
89
91}
92
93template <typename T>
94void
96 const LevelData<T>& a_src,
97 const int a_level,
100 const Interval& a_dstComps,
101 const Interval& a_srcComps,
103 const CopyStrategy& a_fromRegion) const noexcept
104{
105 CH_TIME("AmrMesh::copyData(LD<T>_full)");
106 if (m_verbosity > 5) {
107 pout() << "AmrMesh::copyData(LD<T>_full)" << endl;
108 }
109
110 CH_assert(a_dstComps.size() == a_srcComps.size());
111 CH_assert(a_dst.nComp() > a_dstComps.end());
112 CH_assert(a_src.nComp() > a_srcComps.end());
113
114 if (a_toRealm != a_fromRealm) {
115 const auto id = std::make_pair(a_fromRealm, a_toRealm);
116
117 if (a_toRegion == CopyStrategy::ValidGhost) {
118 CH_assert(a_dst.ghostVect() == m_numGhostCells * IntVect::Unit);
119 }
120 if (a_fromRegion == CopyStrategy::ValidGhost) {
121 CH_assert(a_src.ghostVect() == m_numGhostCells * IntVect::Unit);
122 }
123
125
126 if (a_fromRegion == CopyStrategy::Valid) {
127 if (a_toRegion == CopyStrategy::Valid) {
128 copier = m_validToValidRealmCopiers.at(id)[a_level];
129 }
130 else if (a_toRegion == CopyStrategy::ValidGhost) {
131 copier = m_validToValidGhostRealmCopiers.at(id)[a_level];
132 }
133 else {
134 MayDay::Abort("AmrMesh::copyData - logic bust 1");
135 }
136 }
137 else if (a_fromRegion == CopyStrategy::ValidGhost) {
138 if (a_toRegion == CopyStrategy::Valid) {
139 copier = m_validGhostToValidRealmCopiers.at(id)[a_level];
140 }
141 else if (a_toRegion == CopyStrategy::ValidGhost) {
142 copier = m_validGhostToValidGhostRealmCopiers.at(id)[a_level];
143 }
144 else {
145 MayDay::Abort("AmrMesh::copyData - logic bust 2");
146 }
147 }
148 else {
149 MayDay::Abort("AmrMesh::copyData - logic bust 3");
150 }
151
153 }
154 else {
155 a_src.localCopyTo(a_srcComps, a_dst, a_dstComps);
156 }
157}
158
159template <typename T>
160void
162{
163 CH_TIME("AmrMesh::deallocate(Vector<T*>)");
164 if (m_verbosity > 5) {
165 pout() << "AmrMesh::deallocate(Vector<T*>)" << endl;
166 }
167
168 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
169 delete a_data[lvl];
170 }
171}
172
173template <typename T>
174void
176{
177 CH_TIME("AmrMesh::deallocate(Vector<RefCountedPtr<T>)");
178 if (m_verbosity > 5) {
179 pout() << "AmrMesh::deallocate(Vector<RefCountedPtr<T>)" << endl;
180 }
181
182 for (int lvl = 0; lvl < a_data.size(); lvl++) {
183 // delete &a_data[lvl];
185#if 0
186 if(!a_data[lvl].isNull()){
187 delete &(*a_data[lvl]);
189 }
190
191#endif
192 }
193}
194
195template <typename T>
196void
198{
199 CH_TIME("AmrMesh::deallocate(EBAMRData<T>)");
200 if (m_verbosity > 5) {
201 pout() << "AmrMesh::deallocate(EBAMRData<T>)" << endl;
202 }
203
204 return this->deallocate(a_data.getData());
205}
206
207template <typename T>
208void
210{
211 CH_TIME("AmrMesh::alias(Vector<T*>, Vector<RefCountedPtr<T>)");
212 if (m_verbosity > 5) {
213 pout() << "AmrMesh::alias(Vector<T*>, Vector<RefCountedPtr<T>)" << endl;
214 }
215
216 a_alias.resize(a_data.size());
217
218 for (int lvl = 0; lvl < a_data.size(); lvl++) {
219 a_alias[lvl] = &(*a_data[lvl]);
220 }
221}
222
223template <typename T, typename S>
224void
226{
227 CH_TIME("AmrMesh::alias(Vector<T*>, EBAMRData<S>)");
228 if (m_verbosity > 5) {
229 pout() << "AmrMesh::alias(Vector<T*>, EBAMRData<S>" << endl;
230 }
231
232 return this->alias(a_alias, a_data.getData());
233}
234
235template <typename T>
236void
238{
239 CH_TIME("AmrMesh::allocate(Vector<RefCountedPtr<ParticleData<T> > >, string)");
240 if (m_verbosity > 5) {
241 pout() << "AmrMesh::allocate(Vector<RefCountedPtr<ParticleData<T> > >, string)" << endl;
242 }
243
244 if (!this->queryRealm(a_realm)) {
245 const std::string
246 str = "AmrMesh::allocate(Vector<RefCountedPtr<ParticleData<T> > >, string) - could not find Realm '" + a_realm +
247 "'";
248 MayDay::Abort(str.c_str());
249 }
250
253 "AmrMesh::allocate(Vector<RefCountedPtr<ParticleData<T> > >, string) - only constant box sizes supported for particle methods");
254 }
255
256 a_particles.resize(1 + m_finestLevel);
257
258 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
259 const DisjointBoxLayout& dbl = m_realms[a_realm]->getGrids()[lvl];
260 const ProblemDomain& domain = m_realms[a_realm]->getDomains()[lvl];
261
262 const Real dx = m_dx[lvl];
263
266 }
267}
268
269template <typename T>
270void
272{
273 CH_TIME("AmrMesh::allocate(ParticleContainer<T>, int, string)");
274 if (m_verbosity > 5) {
275 pout() << "AmrMesh::allocate(ParticleContainer<T>, int, string)" << endl;
276 }
277
278 if (!this->queryRealm(a_realm)) {
279 const std::string str = "AmrMesh::allocate(ParticleContainer<T>, int, string) - could not find Realm '" + a_realm +
280 "'";
281 MayDay::Abort(str.c_str());
282 }
283
286 "AmrMesh::allocate(ParticleContainer<T>, int, string) - only constant box sizes are supported for particle methods");
287 }
288
295 m_probLo,
298 a_realm);
299}
300
301template <typename T>
302void
304{
305 CH_TIME("AmrMesh::allocatePointer(Vector<RefCountedPtr<T> >)");
306 if (m_verbosity > 5) {
307 pout() << "AmrMesh::allocatePointer(Vector<RefCountedPtr<T> >)" << endl;
308 }
309
310 a_data.resize(1 + m_finestLevel);
311 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
312 a_data[lvl] = RefCountedPtr<T>(new T());
313 }
314}
315
316template <typename T>
317void
319{
320 CH_TIME("AmrMesh::allocatePointer(Vector<RefCountedPtr<T> >, int)");
321 if (m_verbosity > 5) {
322 pout() << "AmrMesh::allocatePointer(Vector<RefCountedPtr<T> >, int)" << endl;
323 }
324
325 a_data.resize(1 + a_finestLevel);
326
327 for (int lvl = 0; lvl <= a_finestLevel; lvl++) {
328 a_data[lvl] = RefCountedPtr<T>(new T());
329 }
330}
331
332template <typename T>
333void
335{
336 CH_TIME("AmrMesh::allocatePointer(EBAMRData<T>, std::string)");
337 if (m_verbosity > 5) {
338 pout() << "AmrMesh::allocatePointer(EBAMRData<T>, std::string)" << endl;
339 }
340
341 this->allocatePointer(a_data.getData());
342
344}
345
346template <typename T>
347void
349{
350 CH_TIME("AmrMesh::allocatePointer(EBAMRData<T>, std::string, int)");
351 if (m_verbosity > 5) {
352 pout() << "AmrMesh::allocatePointer(EBAMRData<T>, std::string, int)" << endl;
353 }
354
355 this->allocatePointer(a_data.getData(), a_finestLevel);
356
358}
359
360template <class P>
361void
363{
364 CH_TIME("AmrMesh::remapToNewGrids");
365 if (m_verbosity > 5) {
366 pout() << "AmrMesh::remapToNewGrids" << endl;
367 }
368
369 const std::string realm = a_particles.getRealm();
370
371 a_particles.regrid(this->getGrids(realm),
372 this->getDomains(),
373 this->getDx(),
374 this->getRefinementRatios(),
375 this->getValidCells(realm),
376 this->getLevelTiles(realm),
377 a_lmin,
379}
380
381template <class P, class Ret, Ret (P::*MemberFunc)() const>
382void
384 const std::string& a_realm,
385 const phase::which_phase& a_phase,
389 const bool a_forceIrregNGP)
390{
391 CH_TIME("AmrMesh::depositParticles");
392 if (m_verbosity > 5) {
393 pout() << "AmrMesh::depositParticles" << endl;
394 }
395
397
403}
404
405template <class P, typename Ret, Ret (P::*MemberFunc)() const>
406void
408 const std::string& a_realm,
409 const phase::which_phase& a_phase,
410 const ParticleContainer<P>& a_particles) const noexcept
411{
412 CH_TIME("AmrMesh::depositParticles(surface)");
413 if (m_verbosity > 5) {
414 pout() << "AmrMesh::depositParticles(surface)" << endl;
415 }
416
417 CH_assert(a_meshData[0]->nComp() == 1);
418 CH_assert(a_meshData.getRealm() == a_particles.getRealm());
419
420 EBAMRSurfaceDeposition& surfaceDeposition = this->getSurfaceDeposition(a_realm, a_phase);
421
423}
424
425template <class P, class Ret, Ret (P::*MemberFunc)()>
426void
428 const std::string& a_realm,
429 const phase::which_phase& a_phase,
432 const bool a_forceIrregNGP) const
433{
434 CH_TIME("AmrMesh::interpolateParticles(scalar)");
435 if (m_verbosity > 5) {
436 pout() << "AmrMesh::interpolateParticles(scalar)" << endl;
437 }
438
439 CH_assert(a_realm == a_particles.getRealm());
440 CH_assert(a_realm == a_meshScalarField.getRealm());
441
443
445}
446
447template <class P>
448void
450 const phase::which_phase& a_phase,
451 const Real a_tolerance) const
452{
453 CH_TIME("AmrMesh::removeCoveredParticlesIF");
454 if (m_verbosity > 5) {
455 pout() << "AmrMesh::removeCoveredParticlesIF" << endl;
456 }
457
458 CH_assert(!a_particles.isOrganizedByCell());
459
460 // Figure out the implicit function
462
463 switch (a_phase) {
464 case phase::gas: {
465 implicitFunction = m_baseif.at(phase::gas);
466
467 break;
468 }
469 case phase::solid: {
470 implicitFunction = m_baseif.at(phase::solid);
471
472 break;
473 }
474 default: {
475 MayDay::Error("AmrMesh::removeCoveredParticlesIF - logic bust");
476 }
477 }
478
479 // Get the realm where the particles live.
480 const std::string whichRealm = a_particles.getRealm();
481
482 // Go through all particles and remove them if they are less than dx*a_tolerance away from the EB.
483 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
484 const DisjointBoxLayout& dbl = this->getGrids(whichRealm)[lvl];
485 const DataIterator& dit = dbl.dataIterator();
486 const Real dx = this->getDx()[lvl];
487 const Real tol = a_tolerance * dx;
488
489 const int nbox = dit.size();
490#pragma omp parallel for schedule(runtime)
491 for (int mybox = 0; mybox < nbox; mybox++) {
492 const DataIndex& din = dit[mybox];
493
494 List<P>& particles = a_particles[lvl][din].listItems();
495
496 // Check if particles are outside the implicit function.
497 for (ListIterator<P> lit(particles); lit.ok();) {
498 const RealVect& pos = lit().position();
499
500 const Real f = implicitFunction->value(pos);
501
502 if (f > tol) {
503 particles.remove(lit);
504 }
505 else {
506 lit++;
507 }
508 }
509 }
510 }
511}
512
513template <class P>
514void
516 const phase::which_phase& a_phase,
517 const Real a_tolerance) const
518{
519 CH_TIME("AmrMesh::removeCoveredParticlesDiscrete");
520 if (m_verbosity > 5) {
521 pout() << "AmrMesh::removeCoveredParticlesDiscrete" << endl;
522 }
523
524 CH_assert(!a_particles.isOrganizedByCell());
525
526 // Get the realm where the particles live.
527 const std::string whichRealm = a_particles.getRealm();
528
529 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
530 const DisjointBoxLayout& dbl = this->getGrids(whichRealm)[lvl];
531 const DataIterator& dit = dbl.dataIterator();
532 const EBISLayout& ebisl = this->getEBISLayout(whichRealm, a_phase)[lvl];
533 const Real dx = this->getDx()[lvl];
534 const Real tol = a_tolerance * dx;
535
536 const int nbox = dit.size();
537#pragma omp parallel for schedule(runtime)
538 for (int mybox = 0; mybox < nbox; mybox++) {
539 const DataIndex& din = dit[mybox];
540 const EBISBox& ebisBox = ebisl[din];
541 const Box region = ebisBox.getRegion();
542
543 const bool isRegular = ebisBox.isAllRegular();
544 const bool isCovered = ebisBox.isAllCovered();
545 const bool isIrregular = !isRegular && !isCovered;
546
547 List<P>& particles = a_particles[lvl][din].listItems();
548
549 if (isCovered) {
550 particles.clear();
551 }
552 else if (isIrregular) {
553 for (ListIterator<P> lit(particles); lit.ok();) {
554 const P& particle = lit();
555 const RealVect& pos = particle.position();
556
557 // Get the cell where the particle lives.
559
560 CH_assert(region.contains(iv));
561
562 if (ebisBox.isCovered(iv)) {
563 particles.remove(lit);
564 }
565 else if (ebisBox.isIrregular(iv)) {
566
567 // There could potentially be multiple degrees of freedom -- we check if the particle is inside at least one
568 // of the VoFs in the cell.
569
570 bool insideAtLeastOneVoF = false;
571
572 // For each VoF, compute the projection of the particle position onth the EB face. Beacuse the EB normal points outwards, we have a positive
573 // value if the particle lies in the valid region.
574 const std::vector<VolIndex> vofs = ebisBox.getVoFs(iv).stdVector();
575 for (const auto& vof : vofs) {
576 const RealVect ebNormal = ebisBox.normal(vof);
577 const RealVect ebCentroid = m_probLo + Location::position(Location::Cell::Boundary, vof, ebisBox, dx);
578 const Real faceProjection = ebNormal.dotProduct(pos - ebCentroid);
579
580 if (faceProjection >= -tol) {
581 insideAtLeastOneVoF = true;
582
583 break;
584 }
585 }
586
587 if (!insideAtLeastOneVoF) {
588 particles.remove(lit);
589 }
590 else {
591 ++lit;
592 }
593 }
594 else {
595 ++lit;
596 }
597 }
598 }
599 }
600 }
601}
602
603template <class P>
604void
606{
607 CH_TIME("AmrMesh::removeCoveredParticlesVoxels");
608 if (m_verbosity > 5) {
609 pout() << "AmrMesh::removeCoveredParticlesVoxels" << endl;
610 }
611
612 CH_assert(!a_particles.isOrganizedByCell());
613
614 // Get the realm where the particles live.
615 const std::string whichRealm = a_particles.getRealm();
616
617 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
618 const DisjointBoxLayout& dbl = this->getGrids(whichRealm)[lvl];
619 const DataIterator& dit = dbl.dataIterator();
620 const EBISLayout& ebisl = this->getEBISLayout(whichRealm, a_phase)[lvl];
621 const Real dx = this->getDx()[lvl];
622
623 const int nbox = dit.size();
624#pragma omp parallel for schedule(runtime)
625 for (int mybox = 0; mybox < nbox; mybox++) {
626 const DataIndex& din = dit[mybox];
627 const EBISBox& ebisBox = ebisl[din];
628 const Box region = ebisBox.getRegion();
629
630 const bool isRegular = ebisBox.isAllRegular();
631 const bool isCovered = ebisBox.isAllCovered();
632 const bool isIrregular = !isRegular && !isCovered;
633
634 List<P>& particles = a_particles[lvl][din].listItems();
635
636 if (isCovered) {
637 particles.clear();
638 }
639 else if (isIrregular) {
640 for (ListIterator<P> lit(particles); lit.ok();) {
641 const P& particle = lit();
642 const RealVect& pos = particle.position();
643
644 // Get the cell where the particle lives.
645 const RealVect rv = (pos - m_probLo) / dx;
646 const IntVect iv = IntVect(D_DECL(floor(rv[0]), floor(rv[1]), floor(rv[2])));
647
648 CH_assert(region.contains(iv));
649
650 if (ebisBox.isCovered(iv)) {
651 particles.remove(lit);
652 }
653 else {
654 ++lit;
655 }
656 }
657 }
658 }
659 }
660}
661
662template <class P>
663void
666 const phase::which_phase& a_phase,
667 const Real a_tolerance) const
668{
669 CH_TIME("AmrMesh::transferCoveredParticlesIF");
670 if (m_verbosity > 5) {
671 pout() << "AmrMesh::transferCoveredParticlesIF" << endl;
672 }
673
674 CH_assert(!a_particlesFrom.isOrganizedByCell());
675 CH_assert(!a_particlesTo.isOrganizedByCell());
676
677 // Figure out the implicit function
679
680 switch (a_phase) {
681 case phase::gas: {
682 implicitFunction = m_baseif.at(phase::gas);
683
684 break;
685 }
686 case phase::solid: {
687 implicitFunction = m_baseif.at(phase::solid);
688
689 break;
690 }
691 default: {
692 MayDay::Error("AmrMesh::removeCoveredParticlesIF - logic bust");
693
694 break;
695 }
696 }
697
698 // Get the realm where the particles live.
699 const std::string realmFrom = a_particlesFrom.getRealm();
700 const std::string realmTo = a_particlesTo.getRealm();
701
703
704 // Go through all particles and remove them if they are less than dx*a_tolerance away from the EB.
705 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
706 const DisjointBoxLayout& dbl = this->getGrids(realmFrom)[lvl];
707 const DataIterator& dit = dbl.dataIterator();
708 const Real dx = this->getDx()[lvl];
709 const Real tol = a_tolerance * dx;
710
711 const int nbox = dit.size();
712#pragma omp parallel for schedule(runtime)
713 for (int mybox = 0; mybox < nbox; mybox++) {
714 const DataIndex& din = dit[mybox];
715
717 List<P>& particlesTo = a_particlesTo[lvl][din].listItems();
718
719 // Check if particles are outside the implicit function.
720 for (ListIterator<P> lit(particlesFrom); lit.ok();) {
721 const RealVect& pos = lit().position();
722
723 const Real f = implicitFunction->value(pos);
724
725 if (f > tol) {
726 particlesTo.transfer(lit);
727 }
728 else {
729 ++lit;
730 }
731 }
732 }
733 }
734}
735
736template <class P>
737void
740 const phase::which_phase& a_phase,
741 const Real a_tolerance) const
742{
743 CH_TIME("AmrMesh::transferCoveredParticlesDiscrete");
744 if (m_verbosity > 5) {
745 pout() << "AmrMesh::transferCoveredParticlesDiscrete" << endl;
746 }
747
748 CH_assert(!a_particlesFrom.isOrganizedByCell());
749 CH_assert(!a_particlesTo.isOrganizedByCell());
750
751 // Get the realm where the particles live.
752 const std::string realmFrom = a_particlesFrom.getRealm();
753 const std::string realmTo = a_particlesTo.getRealm();
754
756
757 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
758 const DisjointBoxLayout& dbl = this->getGrids(realmFrom)[lvl];
759 const DataIterator& dit = dbl.dataIterator();
760 const EBISLayout& ebisl = this->getEBISLayout(realmFrom, a_phase)[lvl];
761 const Real dx = this->getDx()[lvl];
762 const Real tol = a_tolerance * dx;
763
764 const int nbox = dit.size();
765#pragma omp parallel for schedule(runtime)
766 for (int mybox = 0; mybox < nbox; mybox++) {
767 const DataIndex& din = dit[mybox];
768 const EBISBox& ebisBox = ebisl[din];
769 const Box region = ebisBox.getRegion();
770
771 const bool isRegular = ebisBox.isAllRegular();
772 const bool isCovered = ebisBox.isAllCovered();
773 const bool isIrregular = !isRegular && !isCovered;
774
776 List<P>& particlesTo = a_particlesTo[lvl][din].listItems();
777
778 if (isCovered) {
779 particlesTo.catenate(particlesFrom);
780 }
781 else if (isIrregular) {
782 for (ListIterator<P> lit(particlesFrom); lit.ok();) {
783 const P& particle = lit();
784 const RealVect& pos = particle.position();
785
786 // Get the cell where the particle lives.
788
789 CH_assert(region.contains(iv));
790
791 if (ebisBox.isCovered(iv)) {
792 particlesTo.transfer(lit);
793 }
794 else if (ebisBox.isIrregular(iv)) {
795 // There could potentially be multiple degrees of freedom -- we check if the particle is inside at least one of the VoFs in the cell.
796
797 bool insideAtLeastOneVoF = false;
798
799 // For each VoF, compute the projection of the particle position onth the EB face. Beacuse the EB normal points outwards, we have a positive
800 // value if the particle lies in the valid region.
801 const std::vector<VolIndex> vofs = ebisBox.getVoFs(iv).stdVector();
802 for (const auto& vof : vofs) {
803 const RealVect ebCentroid = m_probLo + Location::position(Location::Cell::Boundary, vof, ebisBox, dx);
804 const RealVect ebNormal = ebisBox.normal(vof);
805 const Real faceProjection = ebNormal.dotProduct(pos - ebCentroid);
806
807 if (faceProjection >= -tol) {
808 insideAtLeastOneVoF = true;
809
810 break;
811 }
812 }
813
814 if (!insideAtLeastOneVoF) {
815 particlesTo.transfer(lit);
816 }
817 else {
818 ++lit;
819 }
820 }
821 else {
822 ++lit;
823 }
824 }
825 }
826 }
827 }
828}
829
830template <class P>
831void
834 const phase::which_phase& a_phase) const
835{
836 CH_TIME("AmrMesh::transferCoveredParticlesVoxels");
837 if (m_verbosity > 5) {
838 pout() << "AmrMesh::transferCoveredParticlesVoxels" << endl;
839 }
840
841 CH_assert(!a_particlesFrom.isOrganizedByCell());
842 CH_assert(!a_particlesTo.isOrganizedByCell());
843
844 // Get the realm where the particles live.
845 const std::string realmFrom = a_particlesFrom.getRealm();
846 const std::string realmTo = a_particlesTo.getRealm();
847
849
850 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
851 const DisjointBoxLayout& dbl = this->getGrids(realmFrom)[lvl];
852 const DataIterator dit = dbl.dataIterator();
853 const EBISLayout& ebisl = this->getEBISLayout(realmFrom, a_phase)[lvl];
854 const Real dx = this->getDx()[lvl];
855
856 const int nbox = dit.size();
857#pragma omp parallel for schedule(runtime)
858 for (int mybox = 0; mybox < nbox; mybox++) {
859 const DataIndex& din = dit[mybox];
860 const EBISBox& ebisBox = ebisl[din];
861 const Box region = ebisBox.getRegion();
862
863 const bool isRegular = ebisBox.isAllRegular();
864 const bool isCovered = ebisBox.isAllCovered();
865 const bool isIrregular = !isRegular && !isCovered;
866
868 List<P>& particlesTo = a_particlesTo[lvl][din].listItems();
869
870 if (isCovered) {
871 particlesTo.catenate(particlesFrom);
872 }
873 else if (isIrregular) {
874 for (ListIterator<P> lit(particlesFrom); lit.ok();) {
875 const P& particle = lit();
876 const RealVect& pos = particle.position();
877
878 // Get the cell where the particle lives.
880
881 CH_assert(region.contains(iv));
882
883 if (ebisBox.isCovered(iv)) {
884 particlesTo.transfer(lit);
885 }
886 else {
887 ++lit;
888 }
889 }
890 }
891 }
892 }
893}
894
895template <class P>
896void
899 const phase::which_phase a_phase,
900 const std::function<void(P&)> a_transferModifier) const noexcept
901{
902 CH_TIME("AmrMesh::transferIrregularParticles");
903 if (m_verbosity > 5) {
904 pout() << "AmrMesh::transferIrregularParticles" << endl;
905 }
906
907 CH_assert(a_dstParticles.getRealm() == a_srcParticles.getRealm());
908
909 const std::string realm = a_dstParticles.getRealm();
910
911 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
912 const DisjointBoxLayout& dbl = this->getGrids(realm)[lvl];
913 const DataIterator& dit = dbl.dataIterator();
914 const EBISLayout& ebisl = this->getEBISLayout(realm, a_phase)[lvl];
915 const ProblemDomain& domain = this->getDomains()[lvl];
916 const RealVect probLo = this->getProbLo();
917 const Real dx = this->getDx()[lvl];
918
919 const int nbox = dit.size();
920#pragma omp parallel for schedule(runtime)
921 for (int mybox = 0; mybox < nbox; mybox++) {
922 const DataIndex& din = dit[mybox];
923
924 const Box cellBox = dbl[din];
925 const EBISBox& ebisbox = ebisl[din];
926
927 List<P>& dstParticles = a_dstParticles[lvl][din].listItems();
928 List<P>& srcParticles = a_srcParticles[lvl][din].listItems();
929
930 for (ListIterator<P> lit(srcParticles); lit.ok();) {
931 P& p = lit();
932
933 const IntVect iv = ParticleOps::getParticleCellIndex(p.position(), probLo, dx);
934
935 if (!(cellBox.contains(iv))) {
936 MayDay::Warning("CD_AmrMeshImplem.H in routine 'transferIrregularParticles' - particle not in box!");
937
938 ++lit;
939 }
940
941 if (ebisbox.isIrregular(iv)) {
942 bool insideAnyVoF = false;
943
944 const Vector<VolIndex> allVoFs = ebisbox.getVoFs(iv);
945
946 for (int i = 0; i < allVoFs.size(); i++) {
947 const VolIndex& vof = allVoFs[i];
948 const RealVect normal = ebisbox.normal(vof);
949 const RealVect ebPos = probLo + Location::position(Location::Cell::Boundary, vof, ebisbox, dx);
950
951 // Note: Can have normal = 0!
952 if ((p.position() - ebPos).dotProduct(normal) >= 0.0) {
953 insideAnyVoF = true;
954 }
955 }
956
957 if (!insideAnyVoF) {
958 const VolIndex& vof = allVoFs[0];
959 const RealVect ebPos = probLo + Location::position(Location::Cell::Boundary, vof, ebisbox, dx);
960
961 p.position() = ebPos;
962
964
965 dstParticles.transfer(lit);
966 }
967 else {
968 ++lit;
969 }
970 }
971 else {
972 ++lit;
973 }
974 }
975 }
976 }
977}
978
979template <class P>
980void
984 const phase::which_phase a_phase,
985 const Real a_tolerance,
986 const bool a_deleteParticles,
987 const std::function<void(P&)> a_nonDeletionModifier) const noexcept
988{
989 CH_TIME("AmrMesh::intersectParticlesRaycastIF");
990 if (m_verbosity > 5) {
991 pout() << "AmrMesh::intersectParticlesRaycastIF" << endl;
992 }
993
994 CH_assert(a_activeParticles.getRealm() == a_ebParticles.getRealm());
995 CH_assert(a_activeParticles.getRealm() == a_domainParticles.getRealm());
996
997 a_ebParticles.clearParticles();
998 a_domainParticles.clearParticles();
999
1000 const std::string whichRealm = a_activeParticles.getRealm();
1001
1002 // Figure out the implicit function
1004
1005 switch (a_phase) {
1006 case phase::gas: {
1007 implicitFunction = m_baseif.at(phase::gas);
1008
1009 break;
1010 }
1011 case phase::solid: {
1012 implicitFunction = m_baseif.at(phase::solid);
1013
1014 break;
1015 }
1016 default: {
1017 MayDay::Error("AmrMesh::intersectParticlesRaycastIF - logic bust");
1018
1019 break;
1020 }
1021 }
1022
1023 // Safety factor to prevent particles falling off the domain if they intersect the high-side of the domain
1024 constexpr Real safety = 1.E-12;
1025
1026 // Level loop -- go through each AMR level
1027 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
1028
1029 // Handle to various grid stuff.
1030 const DisjointBoxLayout& dbl = this->getGrids(whichRealm)[lvl];
1031 const DataIterator& dit = dbl.dataIterator();
1032
1033 const int nbox = dit.size();
1034#pragma omp parallel for schedule(runtime)
1035 for (int mybox = 0; mybox < nbox; mybox++) {
1036 const DataIndex& din = dit[mybox];
1037
1039 List<P>& ebParticles = a_ebParticles[lvl][din].listItems();
1041
1042 for (ListIterator<P> lit(activeParticles); lit.ok();) {
1043 P& particle = lit();
1044
1045 const RealVect newPos = particle.position();
1046 const RealVect oldPos = particle.oldPosition();
1047 const RealVect path = newPos - oldPos;
1048
1049 // Cheap initial tests that allow us to skip some intersections tests.
1050 bool checkEB = false;
1051 bool checkDomain = false;
1052
1053 if (!implicitFunction.isNull()) {
1054 checkEB = true;
1055 }
1056 for (int dir = 0; dir < SpaceDim; dir++) {
1057 const bool outsideLo = newPos[dir] < m_probLo[dir];
1058 const bool outsideHi = newPos[dir] > m_probHi[dir];
1059
1060 if (outsideLo || outsideHi) {
1061 checkDomain = true;
1062 }
1063 }
1064
1065 // Do the intersection tests.
1066 if (checkEB || checkDomain) {
1067
1068 // These are the solution
1069 Real sDomain = std::numeric_limits<Real>::max();
1070 Real sEB = std::numeric_limits<Real>::max();
1071
1072 bool contactDomain = false;
1073 bool contactEB = false;
1074
1075 // Check if the particle intersected the domain. If it did, we compute sDomain such that the intersection
1076 // point with the domain is X = X0 + sDomain * (X1-X0) where X1=newPos and X0=oldPos
1077 if (checkDomain) {
1079 }
1080
1081 // Check if the particle intersected the EB. If it did, we compute sEB such that the intersection
1082 // point with the domain is X = X0 + sEB * (X1-X0) where X1=newPos and X0=oldPos
1083 if (checkEB) {
1085 }
1086
1087 // Particle bumped into something.
1088 if (contactDomain || contactEB) {
1089 if (sEB <= sDomain) { // Crashed with EB "first".
1091
1092 // If we delete the original particles we can just transfer the one we have. Otherwise we create
1093 // a new particle.
1094 if (a_deleteParticles) {
1095 particle.position() = intersectionPos;
1096
1097 ebParticles.transfer(lit);
1098 }
1099 else {
1100 P p = lit();
1101
1102 p.position() = intersectionPos;
1103
1104 ebParticles.add(p);
1105
1107
1108 lit++;
1109 }
1110 }
1111 else {
1112 // Crashed with domain "first". Safety factor is to prevent particles falling off the high side of the domain
1113 const Real sSafety = std::max((Real)0.0, sDomain - safety);
1114
1116
1117 if (a_deleteParticles) {
1118 particle.position() = intersectionPos;
1119
1120 domainParticles.transfer(lit);
1121 }
1122 else {
1123 P p = lit();
1124
1125 p.position() = intersectionPos;
1126
1127 domainParticles.add(p);
1128
1130
1131 ++lit;
1132 }
1133 }
1134 }
1135 else {
1136 ++lit;
1137 }
1138 }
1139 else {
1140 ++lit;
1141 }
1142 }
1143 }
1144 }
1145
1146 // These need to be remapped.
1149}
1150
1151template <class P>
1152void
1156 const phase::which_phase a_phase,
1157 const Real a_bisectionStep,
1158 const bool a_deleteParticles,
1159 const std::function<void(P&)> a_nonDeletionModifier) const noexcept
1160{
1161 CH_TIME("AmrMesh::intersectParticlesBisectIF");
1162 if (m_verbosity > 5) {
1163 pout() << "AmrMesh::intersectParticlesBisectIF" << endl;
1164 }
1165
1166 CH_assert(a_activeParticles.getRealm() == a_ebParticles.getRealm());
1167 CH_assert(a_activeParticles.getRealm() == a_domainParticles.getRealm());
1168
1169 // TLDR: This is pretty much a hard-copy of intersectParticlesRaycastIF, with the exception that the EB intersection
1170 // test is replaced by a bisection test.
1171
1172 a_ebParticles.clearParticles();
1173 a_domainParticles.clearParticles();
1174
1175 const std::string whichRealm = a_activeParticles.getRealm();
1176
1177 // Figure out the implicit function
1179
1180 switch (a_phase) {
1181 case phase::gas: {
1182 implicitFunction = m_baseif.at(phase::gas);
1183
1184 break;
1185 }
1186 case phase::solid: {
1187 implicitFunction = m_baseif.at(phase::solid);
1188
1189 break;
1190 }
1191 default: {
1192 MayDay::Error("AmrMesh::intersectParticlesBisectIF - logic bust");
1193
1194 break;
1195 }
1196 }
1197
1198 // Safety factor to prevent particles falling off the domain if they intersect the high-side of the domain
1199 constexpr Real safety = 1.E-12;
1200
1201 // Level loop -- go through each AMR level
1202 for (int lvl = 0; lvl <= m_finestLevel; lvl++) {
1203
1204 // Handle to various grid stuff.
1205 const DisjointBoxLayout& dbl = this->getGrids(whichRealm)[lvl];
1206 const DataIterator& dit = dbl.dataIterator();
1207
1208 const int nbox = dit.size();
1209#pragma omp parallel for schedule(runtime)
1210 for (int mybox = 0; mybox < nbox; mybox++) {
1211 const DataIndex& din = dit[mybox];
1212
1214 List<P>& ebParticles = a_ebParticles[lvl][din].listItems();
1216
1217 for (ListIterator<P> lit(activeParticles); lit.ok();) {
1218 P& particle = lit();
1219
1220 const RealVect newPos = particle.position();
1221 const RealVect oldPos = particle.oldPosition();
1222 const RealVect path = newPos - oldPos;
1223
1224 // Cheap initial tests that allow us to skip some intersections tests.
1225 bool checkEB = false;
1226 bool checkDomain = false;
1227
1228 if (!implicitFunction.isNull()) {
1229 checkEB = true;
1230 }
1231 for (int dir = 0; dir < SpaceDim; dir++) {
1232 const bool outsideLo = newPos[dir] < m_probLo[dir];
1233 const bool outsideHi = newPos[dir] > m_probHi[dir];
1234
1235 if (outsideLo || outsideHi) {
1236 checkDomain = true;
1237 }
1238 }
1239
1240 // Do the intersection tests.
1241 if (checkEB || checkDomain) {
1242
1243 // These are the solution
1244 Real sDomain = std::numeric_limits<Real>::max();
1245 Real sEB = std::numeric_limits<Real>::max();
1246
1247 bool contactDomain = false;
1248 bool contactEB = false;
1249
1250 // Check if the particle intersected the domain. If it did, we compute sDomain such that the intersection
1251 // point with the domain is X = X0 + sDomain * (X1-X0) where X1=newPos and X0=oldPos
1252 if (checkDomain) {
1254 }
1255
1256 // Check if the particle intersected the EB. If it did, we compute sEB such that the intersection
1257 // point with the domain is X = X0 + sEB * (X1-X0) where X1=newPos and X0=oldPos
1258 if (checkEB) {
1260 }
1261
1262 // Particle bumped into something.
1263 if (contactDomain || contactEB) {
1264 if (sEB <= sDomain) {
1265 // Crashed with EB "first".
1267
1268 if (a_deleteParticles) {
1269 particle.position() = intersectionPos;
1270
1271 ebParticles.transfer(lit);
1272 }
1273 else {
1274 P p = particle;
1275
1276 p.position() = intersectionPos;
1277
1278 ebParticles.add(p);
1279
1280 ++lit;
1281
1283 }
1284 }
1285 else {
1286 // Crashed with domain "first".
1287 const Real sSafety = std::max((Real)0.0, sDomain - safety);
1288
1290
1291 if (a_deleteParticles) {
1292 particle.position() = intersectionPos;
1293
1294 domainParticles.transfer(lit);
1295 }
1296 else {
1297 P p = particle;
1298
1299 p.position() = intersectionPos;
1300
1301 domainParticles.add(p);
1302
1303 ++lit;
1304
1306 }
1307 }
1308 }
1309 else {
1310 ++lit;
1311 }
1312 }
1313 else {
1314 ++lit;
1315 }
1316 }
1317 }
1318 }
1319
1320 // These need to be remapped.
1323}
1324
1325#include <CD_NamespaceFooter.H>
1326
1327#endif
Declaration of core class for handling AMR-related operations (with embedded boundaries)
CoarseFineDeposition
Coarse-fine deposition types (see CD_EBAMRParticleMesh for how these are handled).
Definition CD_CoarseFineDeposition.H:26
CopyStrategy
Enum for distinguishing how we copy data Valid => valid region ValidGhost => valid+ghost region.
Definition CD_CopyStrategy.H:23
DepositionType
Deposition types.
Definition CD_DepositionType.H:23
Declaration of a static class containing some common useful particle routines that would otherwise be...
void removeCoveredParticlesVoxels(ParticleContainer< P > &a_particles, const phase::which_phase &a_phase) const
Function which removes particles from the domain if they fall inside the EB.
Definition CD_AmrMeshImplem.H:605
void transferCoveredParticlesDiscrete(ParticleContainer< P > &a_particlesFrom, ParticleContainer< P > &a_particlesTo, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which transferse particles from one particle container to another if they fall inside the EB...
Definition CD_AmrMeshImplem.H:738
void allocatePointer(Vector< RefCountedPtr< T > > &a_data) const
Allocate pointer but not any memory blocks.
Definition CD_AmrMeshImplem.H:303
const AMRMask & getValidCells(const std::string a_realm) const
Get a map of all valid cells on a specified realm.
Definition CD_AmrMesh.cpp:3225
void copyData(EBAMRData< T > &a_dst, const EBAMRData< T > &a_src, const CopyStrategy &a_toRegion=CopyStrategy::Valid, const CopyStrategy &a_fromRegion=CopyStrategy::Valid) const noexcept
Method for copying from a source container to a destination container. User supplies information abou...
Definition CD_AmrMeshImplem.H:22
void transferCoveredParticlesIF(ParticleContainer< P > &a_particlesFrom, ParticleContainer< P > &a_particlesTo, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which transferse particles from one particle container to another if they fall inside the EB...
Definition CD_AmrMeshImplem.H:664
int m_maxBoxSize
Max box size.
Definition CD_AmrMesh.H:1966
RealVect m_probLo
Domain simulation corner.
Definition CD_AmrMesh.H:1926
void remapToNewGrids(ParticleContainer< P > &a_particles, const int a_lmin, const int a_newFinestLevel) const noexcept
Regrid particle to new grids.
Definition CD_AmrMeshImplem.H:362
void transferIrregularParticles(ParticleContainer< P > &a_dstParticles, ParticleContainer< P > &a_srcParticles, const phase::which_phase a_phase, const std::function< void(P &)> a_transferModifier=[](P &) -> void { return;}) const noexcept
Transfer particles that are on the wrong side of the EB to a different container.
Definition CD_AmrMeshImplem.H:897
bool queryRealm(const std::string a_realm) const
Query if a realm exists.
Definition CD_AmrMesh.cpp:3500
void removeCoveredParticlesIF(ParticleContainer< P > &a_particles, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which removes particles from the domain if they fall inside the EB.
Definition CD_AmrMeshImplem.H:449
int m_verbosity
Verbosity.
Definition CD_AmrMesh.H:1941
void intersectParticlesRaycastIF(ParticleContainer< P > &a_activeParticles, ParticleContainer< P > &a_ebParticles, ParticleContainer< P > &a_domainParticles, const phase::which_phase a_phase, const Real a_tolerance, const bool a_deleteParticles, const std::function< void(P &)> a_nonDeletionModifier=[](P &) -> void { return;}) const noexcept
Particle intersection algorithm based on ray-casting.
Definition CD_AmrMeshImplem.H:981
void removeCoveredParticlesDiscrete(ParticleContainer< P > &a_particles, const phase::which_phase &a_phase, const Real a_tolerance=0.0) const
Function which removes particles from the domain if they fall inside the EB.
Definition CD_AmrMeshImplem.H:515
Vector< Real > m_dx
Level resolutions.
Definition CD_AmrMesh.H:2046
const Vector< DisjointBoxLayout > & getGrids(const std::string a_realm) const
Get the grids.
Definition CD_AmrMesh.cpp:3159
void intersectParticlesBisectIF(ParticleContainer< P > &a_activeParticles, ParticleContainer< P > &a_ebParticles, ParticleContainer< P > &a_domainParticles, const phase::which_phase a_phase, const Real a_bisectionStep, const bool a_deleteParticles, const std::function< void(P &)> a_nonDeletionModifier=[](P &) -> void { return;}) const noexcept
Particle intersection algorithm based on bisection.
Definition CD_AmrMeshImplem.H:1153
int m_finestLevel
Finest level.
Definition CD_AmrMesh.H:1946
const Vector< EBISLayout > & getEBISLayout(const std::string a_realm, const phase::which_phase a_phase) const
Get EBISLayouts for a Realm and phase.
Definition CD_AmrMesh.cpp:3175
void depositParticles(EBAMRCellData &a_meshData, const std::string &a_realm, const phase::which_phase &a_phase, const ParticleContainer< P > &a_particles, const DepositionType a_depositionType, const CoarseFineDeposition a_coarseFineDeposition, const bool a_forceIrregNGP=false)
Deposit scalar particle quantities on the mesh.
Definition CD_AmrMeshImplem.H:383
std::map< phase::which_phase, RefCountedPtr< BaseIF > > m_baseif
Implicit functions.
Definition CD_AmrMesh.H:1845
void allocate(Vector< RefCountedPtr< ParticleData< T > > > &a_particles, const std::string a_realm) const
Template class for generic allocation of particle data.
Definition CD_AmrMeshImplem.H:237
const Vector< Real > & getDx() const
Get spatial resolutions.
Definition CD_AmrMesh.cpp:3104
std::map< std::string, RefCountedPtr< Realm > > m_realms
These are all the Realms.
Definition CD_AmrMesh.H:1840
int m_blockingFactor
Blocking factor.
Definition CD_AmrMesh.H:1981
void interpolateParticles(ParticleContainer< P > &a_particles, const std::string &a_realm, const phase::which_phase &a_phase, const EBAMRCellData &a_meshScalarField, const DepositionType a_interpType, const bool a_forceIrregNGP=false) const
Interpolate mesh data onto a particle position.
Definition CD_AmrMeshImplem.H:427
void transferCoveredParticlesVoxels(ParticleContainer< P > &a_particlesFrom, ParticleContainer< P > &a_particlesTo, const phase::which_phase &a_phase) const
Function which transferse particles from one particle container to another if they fall inside the EB...
Definition CD_AmrMeshImplem.H:832
void alias(Vector< T * > &a_alias, const Vector< RefCountedPtr< T > > &a_data) const
Turn smart-pointer data structure into regular-pointer data structure.
Definition CD_AmrMeshImplem.H:209
EBAMRParticleMesh & getParticleMesh(const std::string a_realm, const phase::which_phase a_phase) const
Get EBAMRParticleMesh operator.
Definition CD_AmrMesh.cpp:3324
void deallocate(Vector< T * > &a_data) const
Deallocate data.
Definition CD_AmrMeshImplem.H:161
const Vector< int > & getRefinementRatios() const
Get refinement ratios.
Definition CD_AmrMesh.cpp:3115
const Vector< RefCountedPtr< LevelTiles > > & getLevelTiles(const std::string a_realm) const
Get the tiled space representation.
Definition CD_AmrMesh.cpp:3241
const Vector< ProblemDomain > & getDomains() const
Get domains.
Definition CD_AmrMesh.cpp:3137
Class for handling particle-mesh operations with AMR.
Definition CD_EBAMRParticleMesh.H:52
class for handling surface deposition of particles with EB and AMR.
Definition CD_EBAMRSurfaceDeposition.H:29
static bool ebIntersectionRaycast(const RefCountedPtr< BaseIF > &a_impFunc, const RealVect &a_oldPos, const RealVect &a_newPos, const Real &a_tolerance, Real &a_s)
Compute the intersection point between a particle path and an implicit function using a ray-casting a...
Definition CD_ParticleOpsImplem.H:221
static bool ebIntersectionBisect(const RefCountedPtr< BaseIF > &a_impFunc, const RealVect &a_oldPos, const RealVect &a_newPos, const Real &a_bisectStep, Real &a_s)
Compute the intersection point between a particle path and an implicit function using a bisection alg...
Definition CD_ParticleOpsImplem.H:177
static IntVect getParticleCellIndex(const RealVect &a_particlePosition, const RealVect &a_probLo, const Real &a_dx) noexcept
Get the cell index corresponding to the particle position.
Definition CD_ParticleOpsImplem.H:30
static bool domainIntersection(const RealVect &a_oldPos, const RealVect &a_newPos, const RealVect &a_probLo, const RealVect &a_probHi, Real &a_s)
Compute the intersection point between a particle path and a domain side.
Definition CD_ParticleOpsImplem.H:131
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
virtual void remap()
Remap particles.
Definition CD_TracerParticleSolverImplem.H:338
virtual void setRealm(const std::string &a_realm)
Set the solver realm.
Definition CD_TracerParticleSolverImplem.H:232
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25
virtual void deposit(EBAMRCellData &a_phi) const noexcept
Deposit particle weight on mesh.
Definition CD_TracerParticleSolverImplem.H:362
virtual void regrid(const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel)
Regrid this solver.
Definition CD_TracerParticleSolverImplem.H:321
RealVect position(const Location::Cell a_location, const VolIndex &a_vof, const EBISBox &a_ebisbox, const Real &a_dx)
Compute the position (ignoring the "origin) of a Vof.
Definition CD_LocationImplem.H:20