Nektar++
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
Nektar::SpatialDomains::BoundaryConditions Class Reference

#include <Conditions.h>

Public Member Functions

 BoundaryConditions (const LibUtilities::SessionReaderSharedPtr &pSession, const MeshGraphSharedPtr &meshGraph)
 
 BoundaryConditions (void)
 
 ~BoundaryConditions (void)
 
const BoundaryRegionCollectionGetBoundaryRegions (void) const
 
void AddBoundaryRegions (const int regionID, BoundaryRegionShPtr &bRegion)
 
const BoundaryConditionCollectionGetBoundaryConditions (void) const
 
void AddBoundaryConditions (const int regionID, BoundaryConditionMapShPtr &bCond)
 
const std::string GetVariable (unsigned int indx)
 
std::map< int, LibUtilities::CommSharedPtrGetBoundaryCommunicators () const
 
const std::map< int, std::string > & GetBoundaryLabels (void) const
 

Protected Attributes

MeshGraphSharedPtr m_meshGraph
 The mesh graph to use for referencing geometry info. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 
BoundaryRegionCollection m_boundaryRegions
 
std::map< int, std::string > m_boundaryLabels
 
BoundaryConditionCollection m_boundaryConditions
 
std::map< int, LibUtilities::CommSharedPtrm_boundaryCommunicators
 

Private Member Functions

void Read (TiXmlElement *conditions)
 Read segments (and general MeshGraph) given TiXmlDocument. More...
 
void ReadBoundaryRegions (TiXmlElement *regions)
 
void ReadBoundaryConditions (TiXmlElement *conditions)
 
void CreateBoundaryComms ()
 

Detailed Description

Definition at line 224 of file Conditions.h.

Constructor & Destructor Documentation

◆ BoundaryConditions() [1/2]

Nektar::SpatialDomains::BoundaryConditions::BoundaryConditions ( const LibUtilities::SessionReaderSharedPtr pSession,
const MeshGraphSharedPtr meshGraph 
)

Constructor - collective on the session's communicator.

Definition at line 48 of file Conditions.cpp.

51 : m_meshGraph(meshGraph), m_session(pSession)
52
53{
54 Read(m_session->GetElement("Nektar/Conditions"));
55}
void Read(TiXmlElement *conditions)
Read segments (and general MeshGraph) given TiXmlDocument.
Definition: Conditions.cpp:220
LibUtilities::SessionReaderSharedPtr m_session
Definition: Conditions.h:273
MeshGraphSharedPtr m_meshGraph
The mesh graph to use for referencing geometry info.
Definition: Conditions.h:272

References m_session, and Read().

◆ BoundaryConditions() [2/2]

Nektar::SpatialDomains::BoundaryConditions::BoundaryConditions ( void  )

Definition at line 57 of file Conditions.cpp.

58{
59}

◆ ~BoundaryConditions()

Nektar::SpatialDomains::BoundaryConditions::~BoundaryConditions ( void  )

Definition at line 61 of file Conditions.cpp.

62{
63}

Member Function Documentation

◆ AddBoundaryConditions()

void Nektar::SpatialDomains::BoundaryConditions::AddBoundaryConditions ( const int  regionID,
BoundaryConditionMapShPtr bCond 
)
inline

Definition at line 249 of file Conditions.h.

251 {
252 m_boundaryConditions[regionID] = bCond;
253 }
BoundaryConditionCollection m_boundaryConditions
Definition: Conditions.h:277

References m_boundaryConditions.

Referenced by Nektar::PulseWaveSystem::SetUpDomainInterfaceBCs().

◆ AddBoundaryRegions()

void Nektar::SpatialDomains::BoundaryConditions::AddBoundaryRegions ( const int  regionID,
BoundaryRegionShPtr bRegion 
)
inline

Definition at line 239 of file Conditions.h.

240 {
241 m_boundaryRegions[regionID] = bRegion;
242 }
BoundaryRegionCollection m_boundaryRegions
Definition: Conditions.h:275

References m_boundaryRegions.

Referenced by Nektar::PulseWaveSystem::SetUpDomainInterfaceBCs().

◆ CreateBoundaryComms()

void Nektar::SpatialDomains::BoundaryConditions::CreateBoundaryComms ( )
private

Create a new communicator for each boundary region. Collective on the session's communicator.

Definition at line 180 of file Conditions.cpp.

181{
182 LibUtilities::CommSharedPtr comm = m_session->GetComm()->GetSpaceComm();
183
184 if (comm->IsSerial())
185 {
186 // Do not try and generate a communicator if we have a serial
187 // communicator. Arises with a FieldConvert communicator when
188 // using --nparts in FieldConvert. Just set communicator to comm
189 // in this case.
190 for (auto &it : m_boundaryRegions)
191 {
192 m_boundaryCommunicators[it.first] = comm;
193 }
194 return;
195 }
196
197 std::set<int> allids = ShareAllBoundaryIDs(m_boundaryRegions, comm);
198
199 for (auto &it : allids)
200 {
201 auto reg_it = m_boundaryRegions.find(it);
202 int this_rank_participates = (reg_it != m_boundaryRegions.end());
203 LibUtilities::CommSharedPtr comm_region =
204 comm->CommCreateIf(this_rank_participates);
205
206 ASSERTL0(bool(comm_region) == bool(this_rank_participates),
207 "Rank should be in communicator but wasn't or is in "
208 "communicator but shouldn't be.");
209
210 if (this_rank_participates)
211 {
212 m_boundaryCommunicators[reg_it->first] = comm_region;
213 }
214 }
215}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
std::map< int, LibUtilities::CommSharedPtr > m_boundaryCommunicators
Definition: Conditions.h:278
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::set< int > ShareAllBoundaryIDs(const BoundaryRegionCollection &boundaryRegions, LibUtilities::CommSharedPtr comm)
Definition: Conditions.cpp:88

References ASSERTL0, m_boundaryCommunicators, m_boundaryRegions, m_session, and Nektar::SpatialDomains::ShareAllBoundaryIDs().

Referenced by Read().

◆ GetBoundaryCommunicators()

std::map< int, LibUtilities::CommSharedPtr > Nektar::SpatialDomains::BoundaryConditions::GetBoundaryCommunicators ( ) const
inline

Definition at line 260 of file Conditions.h.

261 {
263 }

References m_boundaryCommunicators.

Referenced by Nektar::FieldUtils::OutputFileBase::v_Process().

◆ GetBoundaryConditions()

const BoundaryConditionCollection & Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions ( void  ) const
inline

◆ GetBoundaryLabels()

const std::map< int, std::string > & Nektar::SpatialDomains::BoundaryConditions::GetBoundaryLabels ( void  ) const
inline

Definition at line 265 of file Conditions.h.

266 {
267 return m_boundaryLabels;
268 }
std::map< int, std::string > m_boundaryLabels
Definition: Conditions.h:276

References m_boundaryLabels.

Referenced by Nektar::FieldUtils::OutputVtk::OutputFromExpLowOrderMultiBlock().

◆ GetBoundaryRegions()

const BoundaryRegionCollection & Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions ( void  ) const
inline

◆ GetVariable()

const std::string Nektar::SpatialDomains::BoundaryConditions::GetVariable ( unsigned int  indx)
inline

Definition at line 255 of file Conditions.h.

256 {
257 return m_session->GetVariable(indx);
258 }

References m_session.

◆ Read()

void Nektar::SpatialDomains::BoundaryConditions::Read ( TiXmlElement *  conditions)
private

Read segments (and general MeshGraph) given TiXmlDocument.

Collective on the session's communicator.

Definition at line 220 of file Conditions.cpp.

221{
222 ASSERTL0(conditions, "Unable to find CONDITIONS tag in file.");
223
224 TiXmlElement *boundaryRegions =
225 conditions->FirstChildElement("BOUNDARYREGIONS");
226
227 if (boundaryRegions)
228 {
229 ReadBoundaryRegions(conditions);
231 ReadBoundaryConditions(conditions);
232 }
233}
void ReadBoundaryConditions(TiXmlElement *conditions)
Definition: Conditions.cpp:328
void ReadBoundaryRegions(TiXmlElement *regions)
Definition: Conditions.cpp:238

References ASSERTL0, CreateBoundaryComms(), ReadBoundaryConditions(), and ReadBoundaryRegions().

Referenced by BoundaryConditions().

◆ ReadBoundaryConditions()

void Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryConditions ( TiXmlElement *  conditions)
private

Definition at line 328 of file Conditions.cpp.

329{
330 // Protect against multiple reads.
331 if (m_boundaryConditions.size() != 0)
332 {
333 return;
334 }
335
336 // Read REGION tags
337 TiXmlElement *boundaryConditionsElement =
338 conditions->FirstChildElement("BOUNDARYCONDITIONS");
340 boundaryConditionsElement, m_session->GetTimeLevel());
341 ASSERTL0(boundaryConditionsElement,
342 "Boundary conditions must be specified.");
343
344 TiXmlElement *regionElement =
345 boundaryConditionsElement->FirstChildElement("REGION");
346
347 // Read R (Robin), D (Dirichlet), N (Neumann), P (Periodic) C(Cauchy) tags
348 while (regionElement)
349 {
350 BoundaryConditionMapShPtr boundaryConditions =
352
353 int boundaryRegionID;
354 int err = regionElement->QueryIntAttribute("REF", &boundaryRegionID);
355 ASSERTL0(err == TIXML_SUCCESS,
356 "Error reading boundary region reference.");
357
358 ASSERTL0(m_boundaryConditions.count(boundaryRegionID) == 0,
359 "Boundary region '" + std::to_string(boundaryRegionID) +
360 "' appears multiple times.");
361
362 // Find the boundary region corresponding to this ID.
363 std::string boundaryRegionIDStr;
364 std::ostringstream boundaryRegionIDStrm(boundaryRegionIDStr);
365 boundaryRegionIDStrm << boundaryRegionID;
366
367 if (m_boundaryRegions.count(boundaryRegionID) == 0)
368 {
369 regionElement = regionElement->NextSiblingElement("REGION");
370 continue;
371 }
372
373 ASSERTL0(m_boundaryRegions.count(boundaryRegionID) == 1,
374 "Boundary region " +
375 boost::lexical_cast<string>(boundaryRegionID) +
376 " not found");
377
378 // Find the communicator that belongs to this ID
379 LibUtilities::CommSharedPtr boundaryRegionComm =
380 m_boundaryCommunicators[boundaryRegionID];
381
382 TiXmlElement *conditionElement = regionElement->FirstChildElement();
383 std::vector<std::string> vars = m_session->GetVariables();
384
385 while (conditionElement)
386 {
387 // Check type.
388 std::string conditionType = conditionElement->Value();
389 std::string attrData;
390 bool isTimeDependent = false;
391
392 // All have var specified, or else all variables are zero.
393 TiXmlAttribute *attr = conditionElement->FirstAttribute();
394
395 std::vector<std::string>::iterator iter;
396 std::string attrName;
397
398 attrData = conditionElement->Attribute("VAR");
399
400 if (!attrData.empty())
401 {
402 iter = std::find(vars.begin(), vars.end(), attrData);
403 ASSERTL0(
404 iter != vars.end(),
405 (std::string("Cannot find variable: ") + attrData).c_str());
406 }
407
408 if (conditionType == "N")
409 {
410 if (attrData.empty())
411 {
412 // All variables are Neumann and are set to zero.
413 for (auto &varIter : vars)
414 {
415 BoundaryConditionShPtr neumannCondition(
416 MemoryManager<NeumannBoundaryCondition>::
417 AllocateSharedPtr(m_session, "00.0"));
418 (*boundaryConditions)[varIter] = neumannCondition;
419 }
420 }
421 else
422 {
423 if (attr)
424 {
425 std::string equation, userDefined, filename;
426
427 while (attr)
428 {
429
430 attrName = attr->Name();
431
432 if (attrName == "VAR")
433 {
434 // if VAR do nothing
435 }
436 else if (attrName == "USERDEFINEDTYPE")
437 {
438 // Do stuff for the user defined attribute
439 attrData = attr->Value();
440 ASSERTL0(!attrData.empty(),
441 "USERDEFINEDTYPE attribute must have "
442 "associated value.");
443
444 userDefined = attrData;
445 isTimeDependent =
446 boost::iequals(attrData, "TimeDependent");
447 }
448 else if (attrName == "VALUE")
449 {
450 ASSERTL0(attrName == "VALUE",
451 (std::string("Unknown attribute: ") +
452 attrName)
453 .c_str());
454
455 attrData = attr->Value();
456 ASSERTL0(!attrData.empty(),
457 "VALUE attribute must be specified.");
458
459 equation = attrData;
460 }
461 else if (attrName == "FILE")
462 {
463 ASSERTL0(attrName == "FILE",
464 (std::string("Unknown attribute: ") +
465 attrName)
466 .c_str());
467
468 attrData = attr->Value();
469 ASSERTL0(!attrData.empty(),
470 "FILE attribute must be specified.");
471
472 filename = attrData;
473 }
474 else
475 {
476 ASSERTL0(false,
477 (std::string("Unknown boundary "
478 "condition attribute: ") +
479 attrName)
480 .c_str());
481 }
482 attr = attr->Next();
483 }
484
485 BoundaryConditionShPtr neumannCondition(
486 MemoryManager<NeumannBoundaryCondition>::
487 AllocateSharedPtr(m_session, equation,
488 userDefined, filename,
489 boundaryRegionComm));
490 neumannCondition->SetIsTimeDependent(isTimeDependent);
491 (*boundaryConditions)[*iter] = neumannCondition;
492 }
493 else
494 {
495 // This variable's condition is zero.
496 BoundaryConditionShPtr neumannCondition(
497 MemoryManager<NeumannBoundaryCondition>::
498 AllocateSharedPtr(m_session, "0"));
499 (*boundaryConditions)[*iter] = neumannCondition;
500 }
501 }
502 }
503 else if (conditionType == "D")
504 {
505 if (attrData.empty())
506 {
507 // All variables are Dirichlet and are set to zero.
508 for (auto &varIter : vars)
509 {
510 BoundaryConditionShPtr dirichletCondition(
511 MemoryManager<DirichletBoundaryCondition>::
512 AllocateSharedPtr(m_session, "0"));
513 (*boundaryConditions)[varIter] = dirichletCondition;
514 }
515 }
516 else
517 {
518 if (attr)
519 {
520 std::string equation, userDefined, filename;
521
522 while (attr)
523 {
524
525 attrName = attr->Name();
526
527 if (attrName == "VAR")
528 {
529 // if VAR do nothing
530 }
531 else if (attrName == "USERDEFINEDTYPE")
532 {
533
534 // Do stuff for the user defined attribute
535 attrData = attr->Value();
536 ASSERTL0(!attrData.empty(),
537 "USERDEFINEDTYPE attribute must have "
538 "associated value.");
539
540 userDefined = attrData;
541 isTimeDependent =
542 boost::iequals(attrData, "TimeDependent");
543 }
544 else if (attrName == "VALUE")
545 {
546 ASSERTL0(attrName == "VALUE",
547 (std::string("Unknown attribute: ") +
548 attrName)
549 .c_str());
550
551 attrData = attr->Value();
552 ASSERTL0(!attrData.empty(),
553 "VALUE attribute must have associated "
554 "value.");
555
556 equation = attrData;
557
558 if (!boost::iequals(attrData, "0") &&
559 (boost::iequals(userDefined,
560 "WallAdiabatic") ||
561 boost::iequals(userDefined,
562 "WallViscous")))
563 {
564 isTimeDependent = true;
565 }
566 }
567 else if (attrName == "FILE")
568 {
569 ASSERTL0(attrName == "FILE",
570 (std::string("Unknown attribute: ") +
571 attrName)
572 .c_str());
573
574 attrData = attr->Value();
575 ASSERTL0(!attrData.empty(),
576 "FILE attribute must be specified.");
577
578 filename = attrData;
579 }
580 else
581 {
582 ASSERTL0(false,
583 (std::string("Unknown boundary "
584 "condition attribute: ") +
585 attrName)
586 .c_str());
587 }
588 attr = attr->Next();
589 }
590
591 BoundaryConditionShPtr dirichletCondition(
592 MemoryManager<DirichletBoundaryCondition>::
593 AllocateSharedPtr(m_session, equation,
594 userDefined, filename,
595 boundaryRegionComm));
596 dirichletCondition->SetIsTimeDependent(isTimeDependent);
597 (*boundaryConditions)[*iter] = dirichletCondition;
598 }
599 else
600 {
601 // This variable's condition is zero.
602 BoundaryConditionShPtr dirichletCondition(
603 MemoryManager<DirichletBoundaryCondition>::
604 AllocateSharedPtr(m_session, "0"));
605 (*boundaryConditions)[*iter] = dirichletCondition;
606 }
607 }
608 }
609 else if (conditionType == "R") // Read du/dn + PRIMCOEFF u = VALUE
610 {
611 if (attrData.empty())
612 {
613 // All variables are Robin and are set to zero.
614 for (auto &varIter : vars)
615 {
616 BoundaryConditionShPtr robinCondition(
617 MemoryManager<RobinBoundaryCondition>::
618 AllocateSharedPtr(m_session, "0", "0"));
619 (*boundaryConditions)[varIter] = robinCondition;
620 }
621 }
622 else
623 {
624
625 if (attr)
626 {
627 std::string equation1, equation2, userDefined;
628 std::string filename;
629
630 bool primcoeffset = false;
631
632 while (attr)
633 {
634
635 attrName = attr->Name();
636
637 if (attrName == "VAR")
638 {
639 // if VAR do nothing
640 }
641 else if (attrName == "USERDEFINEDTYPE")
642 {
643
644 // Do stuff for the user defined attribute
645 attrData = attr->Value();
646 ASSERTL0(!attrData.empty(),
647 "USERDEFINEDTYPE attribute must have "
648 "associated value.");
649
650 userDefined = attrData;
651 isTimeDependent =
652 boost::iequals(attrData, "TimeDependent");
653 }
654 else if (attrName == "VALUE")
655 {
656
657 attrData = attr->Value();
658 ASSERTL0(!attrData.empty(),
659 "VALUE attributes must have "
660 "associated values.");
661
662 equation1 = attrData;
663 }
664 else if (attrName == "PRIMCOEFF")
665 {
666
667 attrData = attr->Value();
668 ASSERTL0(!attrData.empty(),
669 "PRIMCOEFF attributes must have "
670 "associated values.");
671
672 equation2 = attrData;
673
674 primcoeffset = true;
675 }
676 else if (attrName == "FILE")
677 {
678 attrData = attr->Value();
679 ASSERTL0(!attrData.empty(),
680 "FILE attribute must be specified.");
681
682 filename = attrData;
683 }
684 else
685 {
686 ASSERTL0(false,
687 (std::string("Unknown boundary "
688 "condition attribute: ") +
689 attrName)
690 .c_str());
691 }
692 attr = attr->Next();
693 }
694
695 if (primcoeffset == false)
696 {
697 ASSERTL0(false, "PRIMCOEFF must be specified in a "
698 "Robin boundary condition");
699 }
700 BoundaryConditionShPtr robinCondition(
701 MemoryManager<RobinBoundaryCondition>::
702 AllocateSharedPtr(
703 m_session, equation1, equation2,
704 userDefined, filename, boundaryRegionComm));
705 (*boundaryConditions)[*iter] = robinCondition;
706 }
707 else
708 {
709 // This variable's condition is zero.
710 BoundaryConditionShPtr robinCondition(
711 MemoryManager<RobinBoundaryCondition>::
712 AllocateSharedPtr(m_session, "0", "0"));
713 robinCondition->SetIsTimeDependent(isTimeDependent);
714 (*boundaryConditions)[*iter] = robinCondition;
715 }
716 }
717 }
718 else if (conditionType == "P")
719 {
720 if (attrData.empty())
721 {
722 ASSERTL0(false, "Periodic boundary conditions should "
723 "be explicitely defined");
724 }
725 else
726 {
727
728 if (attr)
729 {
730 std::string userDefined;
731 vector<unsigned int> periodicBndRegionIndex;
732 while (attr)
733 {
734 attrName = attr->Name();
735
736 if (attrName == "VAR")
737 {
738 // if VAR do nothing
739 }
740 else if (attrName == "USERDEFINEDTYPE")
741 {
742 // Do stuff for the user defined attribute
743 attrData = attr->Value();
744 ASSERTL0(!attrData.empty(),
745 "USERDEFINEDTYPE attribute must have "
746 "associated value.");
747
748 userDefined = attrData;
749 isTimeDependent =
750 boost::iequals(attrData, "TimeDependent");
751 }
752 else if (attrName == "VALUE")
753 {
754 attrData = attr->Value();
755 ASSERTL0(!attrData.empty(),
756 "VALUE attribute must have associated "
757 "value.");
758
759 int beg = attrData.find_first_of("[");
760 int end = attrData.find_first_of("]");
761 std::string periodicBndRegionIndexStr =
762 attrData.substr(beg + 1, end - beg - 1);
763 ASSERTL0(
764 beg < end,
765 (std::string("Error reading periodic "
766 "boundary region definition "
767 "for boundary region: ") +
768 boundaryRegionIDStrm.str())
769 .c_str());
770
771 bool parseGood = ParseUtils::GenerateSeqVector(
772 periodicBndRegionIndexStr.c_str(),
773 periodicBndRegionIndex);
774
775 ASSERTL0(
776 parseGood &&
777 (periodicBndRegionIndex.size() == 1),
778 (std::string(
779 "Unable to read periodic boundary "
780 "condition for boundary region: ") +
781 boundaryRegionIDStrm.str())
782 .c_str());
783 }
784 attr = attr->Next();
785 }
786 BoundaryConditionShPtr periodicCondition(
787 MemoryManager<PeriodicBoundaryCondition>::
788 AllocateSharedPtr(periodicBndRegionIndex[0],
789 userDefined,
790 boundaryRegionComm));
791 (*boundaryConditions)[*iter] = periodicCondition;
792 }
793 else
794 {
795 ASSERTL0(false, "Periodic boundary conditions should "
796 "be explicitely defined");
797 }
798 }
799 }
800 else if (conditionType == "C")
801 {
802 ASSERTL0(false, "Cauchy type boundary conditions not "
803 "implemented.");
804 }
805
806 conditionElement = conditionElement->NextSiblingElement();
807 }
808
809 m_boundaryConditions[boundaryRegionID] = boundaryConditions;
810 regionElement = regionElement->NextSiblingElement("REGION");
811 }
812}
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
Definition: ParseUtils.cpp:104
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::shared_ptr< BoundaryConditionMap > BoundaryConditionMapShPtr
Definition: Conditions.h:218
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:447

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::StdRegions::find(), Nektar::ParseUtils::GenerateSeqVector(), Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel(), m_boundaryCommunicators, m_boundaryConditions, m_boundaryRegions, and m_session.

Referenced by Read().

◆ ReadBoundaryRegions()

void Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryRegions ( TiXmlElement *  regions)
private

All elements are of the form: "<B ID="#"> ... </B>", with ? being the element type.

Definition at line 238 of file Conditions.cpp.

239{
240 // ensure boundary regions only read once per class definition
241 if (m_boundaryRegions.size() != 0)
242 {
243 return;
244 }
245
246 TiXmlElement *boundaryRegions =
247 conditions->FirstChildElement("BOUNDARYREGIONS");
249 boundaryRegions, m_session->GetTimeLevel());
250 ASSERTL0(boundaryRegions, "Unable to find BOUNDARYREGIONS block.");
251
252 // See if we have boundary regions defined.
253 TiXmlElement *boundaryRegionsElement =
254 boundaryRegions->FirstChildElement("B");
255
256 while (boundaryRegionsElement)
257 {
258 /// All elements are of the form: "<B ID="#"> ... </B>", with
259 /// ? being the element type.
260 int indx;
261 int err = boundaryRegionsElement->QueryIntAttribute("ID", &indx);
262 ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
263
264 TiXmlNode *boundaryRegionChild = boundaryRegionsElement->FirstChild();
265 // This is primarily to skip comments that may be present.
266 // Comments appear as nodes just like elements.
267 // We are specifically looking for text in the body
268 // of the definition.
269 while (boundaryRegionChild &&
270 boundaryRegionChild->Type() != TiXmlNode::TINYXML_TEXT)
271 {
272 boundaryRegionChild = boundaryRegionChild->NextSibling();
273 }
274
275 ASSERTL0(boundaryRegionChild,
276 "Unable to read variable definition body.");
277 std::string boundaryRegionStr =
278 boundaryRegionChild->ToText()->ValueStr();
279
280 std::string::size_type indxBeg =
281 boundaryRegionStr.find_first_of('[') + 1;
282 std::string::size_type indxEnd =
283 boundaryRegionStr.find_last_of(']') - 1;
284
285 ASSERTL0(indxBeg <= indxEnd,
286 (std::string("Error reading boundary region definition:") +
287 boundaryRegionStr)
288 .c_str());
289
290 std::string indxStr =
291 boundaryRegionStr.substr(indxBeg, indxEnd - indxBeg + 1);
292
293 if (!indxStr.empty())
294 {
295 // Extract the composites from the string and return them in a list.
296 BoundaryRegionShPtr boundaryRegion(
298
299 ASSERTL0(m_boundaryRegions.count(indx) == 0,
300 "Boundary region " + indxStr +
301 " defined more than "
302 "once!");
303
304 m_meshGraph->GetCompositeList(indxStr, *boundaryRegion);
305 if (boundaryRegion->size() > 0)
306 {
307 m_boundaryRegions[indx] = boundaryRegion;
308 }
309
310 // Read optional name as string and save to m_boundaryLabels if
311 // exists
312 std::string name;
313 err = boundaryRegionsElement->QueryStringAttribute("NAME", &name);
314 if (err == TIXML_SUCCESS)
315 {
316 m_boundaryLabels[indx] = name;
317 }
318 }
319
320 boundaryRegionsElement =
321 boundaryRegionsElement->NextSiblingElement("B");
322 }
323}
std::shared_ptr< BoundaryRegion > BoundaryRegionShPtr
Definition: Conditions.h:208

References ASSERTL0, Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel(), m_boundaryLabels, m_boundaryRegions, m_meshGraph, m_session, and CellMLToNektar.pycml::name.

Referenced by Read().

Member Data Documentation

◆ m_boundaryCommunicators

std::map<int, LibUtilities::CommSharedPtr> Nektar::SpatialDomains::BoundaryConditions::m_boundaryCommunicators
protected

◆ m_boundaryConditions

BoundaryConditionCollection Nektar::SpatialDomains::BoundaryConditions::m_boundaryConditions
protected

◆ m_boundaryLabels

std::map<int, std::string> Nektar::SpatialDomains::BoundaryConditions::m_boundaryLabels
protected

Definition at line 276 of file Conditions.h.

Referenced by GetBoundaryLabels(), and ReadBoundaryRegions().

◆ m_boundaryRegions

BoundaryRegionCollection Nektar::SpatialDomains::BoundaryConditions::m_boundaryRegions
protected

◆ m_meshGraph

MeshGraphSharedPtr Nektar::SpatialDomains::BoundaryConditions::m_meshGraph
protected

The mesh graph to use for referencing geometry info.

Definition at line 272 of file Conditions.h.

Referenced by ReadBoundaryRegions().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SpatialDomains::BoundaryConditions::m_session
protected