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

53 : m_meshGraph(meshGraph), m_session(pSession)
54
55{
56 Read(m_session->GetElement("Nektar/Conditions"));
57}
void Read(TiXmlElement *conditions)
Read segments (and general MeshGraph) given TiXmlDocument.
Definition: Conditions.cpp:214
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 59 of file Conditions.cpp.

60{
61}

◆ ~BoundaryConditions()

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

Definition at line 63 of file Conditions.cpp.

64{
65}

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

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

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

215{
216 ASSERTL0(conditions, "Unable to find CONDITIONS tag in file.");
217
218 TiXmlElement *boundaryRegions =
219 conditions->FirstChildElement("BOUNDARYREGIONS");
220
221 if (boundaryRegions)
222 {
223 ReadBoundaryRegions(conditions);
225 ReadBoundaryConditions(conditions);
226 }
227}
void ReadBoundaryConditions(TiXmlElement *conditions)
Definition: Conditions.cpp:322
void ReadBoundaryRegions(TiXmlElement *regions)
Definition: Conditions.cpp:232

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

Referenced by BoundaryConditions().

◆ ReadBoundaryConditions()

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

Definition at line 322 of file Conditions.cpp.

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

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

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