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)=default
 
 ~BoundaryConditions (void)=default
 
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 225 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 46 of file Conditions.cpp.

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

References m_session, and Read().

◆ BoundaryConditions() [2/2]

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

◆ ~BoundaryConditions()

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

Member Function Documentation

◆ AddBoundaryConditions()

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

Definition at line 250 of file Conditions.h.

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

References m_boundaryConditions.

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

◆ AddBoundaryRegions()

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

Definition at line 240 of file Conditions.h.

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

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 170 of file Conditions.cpp.

171{
172 LibUtilities::CommSharedPtr comm = m_session->GetComm()->GetSpaceComm();
173
174 if (comm->IsSerial())
175 {
176 // Do not try and generate a communicator if we have a serial
177 // communicator. Arises with a FieldConvert communicator when
178 // using --nparts in FieldConvert. Just set communicator to comm
179 // in this case.
180 for (auto &it : m_boundaryRegions)
181 {
182 m_boundaryCommunicators[it.first] = comm;
183 }
184 return;
185 }
186
187 std::set<int> allids = ShareAllBoundaryIDs(m_boundaryRegions, comm);
188
189 for (auto &it : allids)
190 {
191 auto reg_it = m_boundaryRegions.find(it);
192 int this_rank_participates = (reg_it != m_boundaryRegions.end());
193 LibUtilities::CommSharedPtr comm_region =
194 comm->CommCreateIf(this_rank_participates);
195
196 ASSERTL0(bool(comm_region) == bool(this_rank_participates),
197 "Rank should be in communicator but wasn't or is in "
198 "communicator but shouldn't be.");
199
200 if (this_rank_participates)
201 {
202 m_boundaryCommunicators[reg_it->first] = comm_region;
203 }
204 }
205}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
std::map< int, LibUtilities::CommSharedPtr > m_boundaryCommunicators
Definition: Conditions.h:279
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:78

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 261 of file Conditions.h.

262 {
264 }

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 266 of file Conditions.h.

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

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 256 of file Conditions.h.

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

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 210 of file Conditions.cpp.

211{
212 ASSERTL0(conditions, "Unable to find CONDITIONS tag in file.");
213
214 TiXmlElement *boundaryRegions =
215 conditions->FirstChildElement("BOUNDARYREGIONS");
216
217 if (boundaryRegions)
218 {
219 ReadBoundaryRegions(conditions);
221 ReadBoundaryConditions(conditions);
222 }
223}
void ReadBoundaryConditions(TiXmlElement *conditions)
Definition: Conditions.cpp:318
void ReadBoundaryRegions(TiXmlElement *regions)
Definition: Conditions.cpp:228

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

Referenced by BoundaryConditions().

◆ ReadBoundaryConditions()

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

Definition at line 318 of file Conditions.cpp.

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

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 228 of file Conditions.cpp.

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

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 277 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 273 of file Conditions.h.

Referenced by ReadBoundaryRegions().

◆ m_session

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