Nektar++
PulseWaveSystem.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PulseWaveSystem.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Generic timestepping for Pulse Wave Solver
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#include <iostream>
37
43using namespace std;
44
45namespace Nektar
46{
47
48/**
49 * @class PulseWaveSystem
50 *
51 * Initialises the arterial subdomains in m_vessels and sets up
52 * all domain-linking conditions (bifurcations, intefaces/junctions,
53 * merging flows). Detects the network structure and assigns
54 * boundary conditons. Also provides the underlying timestepping
55 * framework for pulse wave solvers including the general
56 * timestepping routines.
57 */
58
59/**
60 * Processes SolverInfo parameters from the session file and sets
61 * up timestepping-specific code.
62 *
63 * @param m_Session Session object to read parameters from.
64 */
68 : UnsteadySystem(pSession, pGraph)
69{
70}
71
72/**
73 * Destructor
74 */
76{
77}
78
79/**
80 * Initialisation routine for multidomain solver. Sets up the
81 * expansions for every arterial segment (m_vessels) and for one
82 * complete field m_outfield which is needed to write the
83 * postprocessing output. Also determines which upwind strategy
84 * is used (currently only upwinding scheme available) and reads
85 * blodd flow specific parameters from the inputfile
86 *
87 */
88void PulseWaveSystem::v_InitObject(bool DeclareField)
89{
90 // Initialise base class
91 UnsteadySystem::v_InitObject(DeclareField);
92
93 // Read the geometry and the expansion information
94 m_domain = m_graph->GetDomain();
95 m_nDomains = m_domain.size();
96
97 // Determine projectiontype
98 ASSERTL0(m_session->MatchSolverInfo("Projection", "DisContinuous"),
99 "Pulse solver only set up for Discontinuous projections");
101 ASSERTL0(
102 m_graph->GetMeshDimension() == 1,
103 "Pulse wave solver only set up for expansion dimension equal to 1");
104
105 size_t i;
106 m_nVariables = m_session->GetVariables().size();
107
109 m_vessels =
111
112 /* If the PressureArea property is not specified, default to the Beta
113 * law; it's the most well known and this way previous examples that did
114 * not specify the tube law will still work.
115 */
116 if (m_session->DefinesSolverInfo("PressureArea"))
117 {
119 m_session->GetSolverInfo("PressureArea"), m_vessels, m_session);
120 }
121 else
122 {
124 "Beta", m_vessels, m_session);
125 }
126
128
132
134
135 // Set up domains and put geometry to be only one space dimension.
136 size_t cnt = 0;
137 bool SetToOneSpaceDimension = true;
138
139 if (m_session->DefinesCmdLineArgument("SetToOneSpaceDimension"))
140 {
141 std::string cmdline = m_session->GetCmdLineArgument<std::string>(
142 "SetToOneSpaceDimension");
143 if (boost::to_upper_copy(cmdline) == "FALSE")
144 {
145 SetToOneSpaceDimension = false;
146 }
147 }
148
149 // get parallel communicators as required
150 map<int, LibUtilities::CommSharedPtr> domComm;
151 GetCommArray(domComm);
152
153 SetUpDomainInterfaceBCs(Allbcs, domComm);
154
155 for (auto &d : m_domOrder)
156 {
157 for (size_t j = 0; j < m_nVariables; ++j)
158 {
159 m_vessels[cnt++] =
161 m_session, m_graph, m_domain[d], Allbcs,
162 m_session->GetVariable(j), domComm[d],
163 SetToOneSpaceDimension);
164 }
165 }
166
167 // Reset coeff and phys space to be continuous over all domains
168 int totcoeffs = 0;
169 int totphys = 0;
171 for (i = 0; i < m_nDomains; ++i)
172 {
173 totcoeffs += m_vessels[i * m_nVariables]->GetNcoeffs();
174
175 m_fieldPhysOffset[i] = totphys;
176 totphys += m_vessels[i * m_nVariables]->GetTotPoints();
177 }
178 m_fieldPhysOffset[m_nDomains] = totphys;
179
180 for (size_t n = 0; n < m_nVariables; ++n)
181 {
182 Array<OneD, NekDouble> coeffs(totcoeffs, 0.0);
183 Array<OneD, NekDouble> phys(totphys, 0.0);
184 Array<OneD, NekDouble> tmpcoeffs, tmpphys;
185
186 m_vessels[n]->SetCoeffsArray(coeffs);
187 m_vessels[n]->SetPhysArray(phys);
188
189 int cnt = m_vessels[n]->GetNcoeffs();
190 int cnt1 = m_vessels[n]->GetTotPoints();
191
192 for (i = 1; i < m_nDomains; ++i)
193 {
194 m_vessels[i * m_nVariables + n]->SetCoeffsArray(tmpcoeffs =
195 coeffs + cnt);
196 m_vessels[i * m_nVariables + n]->SetPhysArray(tmpphys =
197 phys + cnt1);
198 cnt += m_vessels[i * m_nVariables + n]->GetNcoeffs();
199 cnt1 += m_vessels[i * m_nVariables + n]->GetTotPoints();
200 }
201 }
202
203 // If Discontinuous Galerkin determine upwinding method to use
204 for (size_t i = 0; i < SIZE_UpwindTypePulse; ++i)
205 {
206 bool match;
207 m_session->MatchSolverInfo("UPWINDTYPEPULSE", UpwindTypeMapPulse[i],
208 match, false);
209 if (match)
210 {
212 break;
213 }
214 }
215
216 // Load blood density and external pressure
217 m_session->LoadParameter("rho", m_rho, 0.5);
218 m_session->LoadParameter("pext", m_pext, 0.0);
219
220 int nq = 0;
221 /*
222 * Gets the Material Properties of each arterial segment
223 * specified in the inputfile from section MaterialProperties
224 * Also gets the Area at static equilibrium A_0 specified in the
225 * inputfile.
226 *
227 * Having found these points also extract the values at the
228 * trace points and the normal direction consistent with the
229 * left adjacent definition of Fwd and Bwd
230 */
239
240 int omega = 0;
241 for (auto &d : m_domOrder)
242 {
243 nq = m_vessels[2 * omega]->GetNpoints();
244 m_fields[0] = m_vessels[2 * omega];
245
246 m_beta[omega] = Array<OneD, NekDouble>(nq);
247 GetFunction("MaterialProperties")
248 ->Evaluate("beta", m_beta[omega], m_time, d);
249
250 // If the input file doesn't specify viscoelasticity, set it to zero.
251 m_gamma[omega] = Array<OneD, NekDouble>(nq);
252
253 if (m_session->DefinesFunction("Viscoelasticity"))
254 {
255 GetFunction("Viscoelasticity")
256 ->Evaluate("gamma", m_gamma[omega], m_time, d);
257 }
258 else
259 {
260 for (int j = 0; j < nq; ++j)
261 {
262 m_gamma[omega][j] = 0;
263 }
264 }
265
266 /* If the input file doesn't specify strain-stiffening, set it to
267 * 0.5 for elastic behaviour.
268 */
269 m_alpha[omega] = Array<OneD, NekDouble>(nq);
270
271 if (m_session->DefinesFunction("StrainStiffening"))
272 {
273 GetFunction("StrainStiffening")
274 ->Evaluate("alpha", m_alpha[omega], m_time, d);
275 }
276 else
277 {
278 for (int j = 0; j < nq; ++j)
279 {
280 m_alpha[omega][j] = 0.5;
281 }
282 }
283
284 m_A_0[omega] = Array<OneD, NekDouble>(nq);
285 GetFunction("A_0")->Evaluate("A_0", m_A_0[omega], m_time, d);
286
287 int nqTrace = GetTraceTotPoints();
288
289 m_beta_trace[omega] = Array<OneD, NekDouble>(nqTrace);
290 m_fields[0]->ExtractTracePhys(m_beta[omega], m_beta_trace[omega]);
291
292 m_A_0_trace[omega] = Array<OneD, NekDouble>(nqTrace);
293 m_fields[0]->ExtractTracePhys(m_A_0[omega], m_A_0_trace[omega]);
294
295 m_alpha_trace[omega] = Array<OneD, NekDouble>(nqTrace);
296 m_fields[0]->ExtractTracePhys(m_alpha[omega], m_alpha_trace[omega]);
297
298 if (SetToOneSpaceDimension)
299 {
300 m_trace_fwd_normal[omega] = Array<OneD, NekDouble>(nqTrace, 0.0);
301
302 MultiRegions::ExpListSharedPtr trace = m_fields[0]->GetTrace();
303 int nelmt_trace = trace->GetExpSize();
304
305 Array<OneD, Array<OneD, NekDouble>> normals(nelmt_trace);
306
307 for (int i = 0; i < nelmt_trace; ++i)
308 {
309 normals[i] = m_trace_fwd_normal[omega] + i;
310 }
311
312 // need to set to 1 for consistency since boundary
313 // conditions may not have coordim = 1
314 trace->GetExp(0)->GetGeom()->SetCoordim(1);
315
316 trace->GetNormals(normals);
317 }
318 omega++;
319 }
320
322}
323
325 map<int, LibUtilities::CommSharedPtr> &retval)
326{
327
328 // set up defualt serial communicator
329 std::string def("default");
330 char *argv = new char[def.length() + 1];
331 std::strcpy(argv, def.c_str());
332 LibUtilities::CommSharedPtr serialComm =
334
335 size_t nprocs = m_comm->GetSize();
336
337 if (nprocs == 1) // serial case
338 {
339 for (auto &d : m_domain)
340 {
341 // serial communicator
342 retval[d.first] = serialComm;
343 // fill out m_domOrder set for serial case
344 m_domOrder.push_back(d.first);
345 }
346 }
347 else // parallel case
348 {
349 size_t rank = m_comm->GetRank();
350
351 // Fill array with domain details
352 size_t dmax = 0;
353 for (auto &d : m_domain)
354 {
355 dmax = ((size_t)d.first > dmax) ? d.first : dmax;
356 }
357 dmax += 1;
358
359 // get domain maximum id + 1;
360 m_comm->AllReduce(dmax, LibUtilities::ReduceMax);
361
362 Array<OneD, int> commtmp(dmax, 0);
363 for (auto &d : m_domain)
364 {
365 commtmp[d.first] = 1;
366 }
367
368 // identfy which domains are split between two procs
369 m_comm->AllReduce(commtmp, LibUtilities::ReduceSum);
370
371 Array<OneD, bool> DoneCreate(nprocs, false);
372
373 // setup communicators for domain partitions greater than 2
374 for (size_t d = 0; d < dmax; ++d)
375 {
376 // set up communicator for cases where more than 2 procs
377 // share domain
378 if (commtmp[d] > 2)
379 {
381 int flag = 0;
382 // set flag to 10 if d is on this processors domain
383 for (auto &d1 : m_domain)
384 {
385 if ((size_t)d1.first == d)
386 {
387 flag = 10;
388 }
389 }
390
391 newcomm = m_comm->CommCreateIf(flag);
392
393 // set up communicator if domain on this processor
394 for (auto &d1 : m_domain)
395 {
396 if ((size_t)d1.first == d)
397 {
398 retval[d] = newcomm;
399 }
400 }
401
402 // turn off commtmp flag since should now be setup
403 commtmp[d] = -1;
404 }
405 }
406
407 // we now should have a list of the number of each
408 // domains/vessel that is split over two processors
409
410 // set up serial communicators for domains on only one
411 // processor and set up map of domains that require parallel
412 // communicator and number of processors in parallel
413
414 // communicator
415 map<int, int> SharedProc;
416 Array<OneD, int> commtmp1(dmax, 0);
417
418 for (auto &d : m_domain)
419 {
420 if (commtmp[d.first] == 1) // set to serial communicator
421 {
422 retval[d.first] = serialComm;
423 }
424 else if (commtmp[d.first] == 2) // shared by 2 procs
425 {
426 SharedProc[d.first] = 1;
427 commtmp1[d.first] = rank;
428 }
429 else // do nothing
430 {
431 }
432 }
433
434 // find the maximum number of communicator on all processors
435 size_t nShrProc = SharedProc.size();
436 int commMax = nShrProc;
437 m_comm->AllReduce(commMax, LibUtilities::ReduceMax);
438
439 // find out processor id of each communicator by adding rank
440 // and then subtracting local rank
441
442 m_comm->AllReduce(commtmp1, LibUtilities::ReduceSum);
443
444 // reset SharedProc to hold shared proc id.
445 for (auto &d : SharedProc)
446 {
447 SharedProc[d.first] = commtmp1[d.first] - rank;
448 }
449
450 // now work out a comm ordering;
451 Array<OneD, int> procs(2);
452 int commcall = 0;
453 map<int, pair<int, int>> CreateComm;
454
455 bool search = true;
456 Array<OneD, int> setdone(3);
457
458 while (search)
459 {
460 size_t dorank = 0; // search index over processors
461 set<int>
462 proclist; // has proc been identified for a comm at this level
463 int flag = 100; // colour for communicators in this search level
464
465 while (dorank < nprocs)
466 {
467 int sharedproc = -1;
468 int ncomms = 0;
469
470 if (DoneCreate[dorank] == false)
471 {
472 // check to see if proc even has a shared domain
473 ncomms = (rank == dorank) ? nShrProc : 0;
474 m_comm->AllReduce(ncomms, LibUtilities::ReduceMax);
475
476 // set DoneCreate to true if n comms required on
477 // proc
478 DoneCreate[dorank] = (ncomms == 0) ? true : false;
479 }
480
481 // requires communicator and needs to create one.
482 if (ncomms)
483 {
484 bool doneCreateOnDoRank = false;
485 bool doneCreateOnSharedProc = false;
486 bool createdComm = false;
487
488 // check not already callng communicator on
489 // this processor for this search sweep
490 if ((proclist.count(dorank) == 0))
491 {
492
493 // give all processors rank id of shared proc
494 if (rank == dorank)
495 {
496 for (auto &d : SharedProc)
497 {
498 if (CreateComm.count(d.first) == 0)
499 {
500 sharedproc = SharedProc[d.first];
501 break;
502 }
503 }
504
505 ASSERTL1(sharedproc != -1,
506 "Failed to fine a new "
507 "processor to set up another "
508 "communicator");
509
510 if (proclist.count(sharedproc) == 0)
511 {
512 // save all communicators on this for
513 // split comm setup
514 for (auto &d : SharedProc)
515 {
516 if (d.second == sharedproc)
517 {
518 CreateComm[d.first] = pair<int, int>(
519 commcall, flag + d.first);
520 createdComm = true;
521 }
522 }
523 }
524
525 // set bool to mark as done on dorank
526 if (CreateComm.size() == nShrProc)
527 {
528 doneCreateOnDoRank = true;
529 }
530 }
531 else
532 {
533 sharedproc = -1;
534 }
535
536 m_comm->AllReduce(sharedproc, LibUtilities::ReduceMax);
537 ASSERTL1(
538 sharedproc < nprocs,
539 "Error in obtaining rank "
540 "shared by domain. Sharedproc value is greater "
541 "than number of procs");
542
543 if (proclist.count(sharedproc) == 0)
544 {
545 if ((int)rank == sharedproc)
546 {
547 // save all communicators on this for
548 // split comm setup
549 for (auto &d : SharedProc)
550 {
551 if (d.second == (int)dorank)
552 {
553 CreateComm[d.first] = pair<int, int>(
554 commcall, flag + d.first);
555 createdComm = true;
556 }
557 }
558
559 if (CreateComm.size() == nShrProc)
560 {
561 doneCreateOnSharedProc = true;
562 }
563 }
564 }
565 }
566
567 // mark DoneCreate as if complete
568 setdone[0] = (doneCreateOnDoRank) ? 1 : 0;
569 setdone[1] = (doneCreateOnSharedProc) ? 1 : 0;
570 setdone[2] = (createdComm) ? 1 : 0;
571 m_comm->AllReduce(setdone, LibUtilities::ReduceMax);
572 DoneCreate[dorank] = (bool)setdone[0];
573 if (sharedproc != -1)
574 {
575 DoneCreate[sharedproc] = (bool)setdone[1];
576 }
577 if (setdone[2]) // created new communicator
578 {
579 proclist.insert(dorank);
580 proclist.insert(sharedproc);
581 }
582 }
583
584 dorank++;
585 }
586
587 // have we found all comms on all processors.
588 size_t i;
589 for (i = 0; i < nprocs; ++i)
590 {
591 if (DoneCreate[i] == false)
592 {
593 break;
594 }
595 }
596
597 if (i == nprocs)
598 {
599 search = false;
600 }
601 else
602 {
603 commcall++;
604 ASSERTL0(commcall < 4 * commMax,
605 "Failed to find sub communicators "
606 "pattern in 4*commMax searches");
607 }
608 }
609
610 ASSERTL1(CreateComm.size() == nShrProc,
611 "Have not created communicators for all shared procs");
612
613 // determine maxmimum number of CreateComm size
614 size_t maxCreateComm = CreateComm.size();
615 m_comm->AllReduce(maxCreateComm, LibUtilities::ReduceMax);
616
617 // loop over CreateComm list
618 for (size_t i = 0; i < maxCreateComm; ++i)
619 {
621 int flag = 0;
622
623 for (auto &d : CreateComm)
624 {
625 if (d.second.first == (int)i)
626 {
627 flag = d.second.second;
628 }
629 }
630
631 newcomm = m_comm->CommCreateIf(flag);
632
633 for (auto &d : CreateComm)
634 {
635 if (d.second.first == (int)i)
636 {
637 retval[d.first] = newcomm;
638 }
639 }
640 }
641
642 // Finally need to re-order domain list so that we call the
643 // communicators in a non-locking manner. This is done by
644 // uniquely numbering each shared proc/comm interface and
645 // ensuring the domains are ordered so that we call this
646 // numbering in an increasing order
647
648 Array<OneD, NekDouble> numShared(nprocs, 0.0);
649 // determine the number of communicators on this rank where
650 // the share proc has a higher rank id
651 numShared[rank] = nShrProc;
652 m_comm->AllReduce(numShared, LibUtilities::ReduceSum);
653 int totShared = (int)Vmath::Vsum(nprocs, numShared, 1);
654
655 if (totShared)
656 {
657 size_t cnt = 0;
658 Array<OneD, NekDouble> numOffset(nprocs, 0.0);
659 for (auto &s : SharedProc)
660 {
661 if (s.second > (int)rank)
662 {
663 cnt++;
664 }
665 }
666
667 numOffset[rank] = cnt;
668 m_comm->AllReduce(numOffset, LibUtilities::ReduceMax);
669
670 // make numShared into a cumulative list
671 for (size_t i = 1; i < nprocs; ++i)
672 {
673 numShared[i] += numShared[i - 1];
674 numOffset[i] += numOffset[i - 1];
675 }
676 for (size_t i = nprocs - 1; i > 0; --i)
677 {
678 numShared[i] = numShared[i - 1];
679 numOffset[i] = numOffset[i - 1];
680 }
681 numShared[0] = 0;
682 numOffset[0] = 0;
683
684 Array<OneD, NekDouble> shareddom(totShared, -1.0);
685 cnt = 0; // collect a list of domains that each shared commm is
686 // attached to
687 for (auto &s : SharedProc)
688 {
689 shareddom[numShared[rank] + cnt] = s.first;
690 ++cnt;
691 }
692 m_comm->AllReduce(shareddom, LibUtilities::ReduceMax);
693
694 // define a numbering scheme
695 Array<OneD, NekDouble> sharedid(totShared, -1.0);
696 cnt = 0;
697 size_t cnt1 = 0;
698 for (auto &s : SharedProc)
699 {
700 if (s.second > (int)rank)
701 {
702 sharedid[numShared[rank] + cnt] = numOffset[rank] + cnt1;
703
704 // find shared proc offset by matching the domain ids.
705 size_t j;
706 for (j = 0; j < maxCreateComm; ++j)
707 {
708 if ((numShared[s.second] + j < totShared) &&
709 (shareddom[numShared[s.second] + j] == s.first))
710 {
711 break;
712 }
713 }
714
715 sharedid[numShared[s.second] + j] = numOffset[rank] + cnt1;
716 cnt1++;
717 }
718 cnt++;
719 }
720
721 // now communicate ids to all procesors.
722 m_comm->AllReduce(sharedid, LibUtilities::ReduceMax);
723
724 if (rank == 0)
725 {
726 for (size_t i = 0; i < (size_t)totShared; ++i)
727 {
728 ASSERTL1(sharedid[i] != -1.0,
729 "Failed to number shared proc uniquely");
730 }
731 }
732
733 cnt = 0;
734 int maxoffset =
735 (int)Vmath::Vmax(nShrProc, &sharedid[numShared[rank]], 1);
736 m_comm->AllReduce(maxoffset, LibUtilities::ReduceMax);
737 maxoffset++;
738
739 // make a map relating the order of the SharedProc to the domain id
740 map<int, int> ShrToDom;
741 cnt = 0;
742 for (auto &s : SharedProc)
743 {
744 ShrToDom[cnt] = s.first;
745 ++cnt;
746 }
747
748 // Set up a set the domain ids listed in the ordering of
749 // the SharedProc unique numbering
750 set<int> doneDom;
751 NekDouble maxdom = Vmath::Vmax(totShared, shareddom, 1);
752 maxdom++;
753 for (size_t i = 0; i < nShrProc; ++i)
754 {
755 int minId =
756 Vmath::Imin(nShrProc, &sharedid[numShared[rank]], 1);
757 sharedid[numShared[rank] + minId] += maxoffset;
758
759 m_domOrder.push_back(ShrToDom[minId]);
760 doneDom.insert(ShrToDom[minId]);
761 }
762
763 // add all the other
764 for (auto &d : m_domain)
765 {
766 if (doneDom.count(d.first) == 0)
767 {
768 m_domOrder.push_back(d.first);
769 }
770 }
771 }
772 }
773}
774
776{
777 map<int, std::vector<InterfacePointShPtr>> VidToDomain;
778
779 // First set up domain to determine how many vertices meet at a
780 // vertex. This allows us to
781
782 /* Loop over domain and find out if we have any undefined
783 * boundary conditions representing interfaces. If so make a
784 * map based around vid and storing the domains that are
785 * part of interfaces.
786 */
787 for (size_t omega = 0; omega < m_nDomains; ++omega)
788 {
789 size_t vesselID = omega * m_nVariables;
790
791 for (size_t i = 0; i < (m_vessels[vesselID]->GetBndConditions()).size();
792 ++i)
793 {
794 if (m_vessels[vesselID]->GetBndConditions()[i]->GetUserDefined() ==
795 "Interface")
796 {
797 // Get Vid of interface
798 int vid = m_vessels[vesselID]
799 ->UpdateBndCondExpansion(i)
800 ->GetExp(0)
801 ->GetGeom()
802 ->GetVid(0);
803
805 m_vessels[vesselID]->GetTrace();
807
808 bool finish = false;
809 // find which elmt, the local vertex and the data
810 // offset of point
811 size_t nExp = m_vessels[vesselID]->GetExpSize();
812 for (size_t n = 0; n < nExp; ++n)
813 {
814 for (size_t p = 0; p < 2; ++p)
815 {
816
817 if (m_vessels[vesselID]
818 ->GetTraceMap()
819 ->GetElmtToTrace()[n][p]
821 ->GetGeom()
822 ->GetVid(0) == vid)
823 {
824 int eid = m_vessels[vesselID]
825 ->GetTraceMap()
826 ->GetElmtToTrace()[n][p]
827 ->GetElmtId();
828
829 int tid = m_vessels[vesselID]
830 ->GetTrace()
831 ->GetCoeff_Offset(eid);
832
833 Ipt = MemoryManager<
834 InterfacePoint>::AllocateSharedPtr(vid, omega,
835 n, p, tid,
836 i);
837
838 finish = true;
839 break;
840 }
841 }
842 if (finish == true)
843 {
844 break;
845 }
846 }
847
848 VidToDomain[vid].push_back(Ipt);
849 }
850 }
851 }
852
853 // Set up comms for parallel cases
854 map<int, int> domId;
855
856 size_t cnt = 0;
857 for (auto &d : m_domain)
858 {
859 domId[cnt] = d.first;
860 ++cnt;
861 }
862 ASSERTL1(domId.size() == m_nDomains, "Number of domains do not match");
863
864 Gs::gs_data *gscomm;
865 int nvid2dom = VidToDomain.size();
866 Array<OneD, long> tmp(nvid2dom);
867 Array<OneD, NekDouble> nvid(nvid2dom, 0.0);
868
869 cnt = 0;
870 for (auto &v : VidToDomain)
871 {
872 tmp[cnt] = v.first + 1;
873 for (size_t i = 0; i < v.second.size(); ++i)
874 {
875 nvid[cnt] += 1.0;
876 }
877 cnt++;
878 }
879
880 const bool verbose = m_session->DefinesCmdLineArgument("verbose");
881
882 gscomm = Gs::Init(tmp, m_comm->GetRowComm(), verbose);
883 Gs::Gather(nvid, Gs::gs_add, gscomm);
884
885 // Loop over domains and identify how many vessels at a
886 // bifurcation use the beginning of the element at the
887 // junction. This distinguishes a bifurcation and merging junction
888 // since we implicitly assume elements are ordered left ot right
889 // and so at a bifurcation we have two elements meeting this point
890 // at the first vertex
891 Array<OneD, NekDouble> nbeg(nvid2dom, 0.0);
892 cnt = 0;
893 for (auto &v : VidToDomain)
894 {
895 if (nvid[cnt] == 3.0)
896 {
897 for (size_t i = 0; i < v.second.size(); ++i)
898 {
899 // for bifurcations and merging junctions store how
900 // many points at interface are at the beginning or
901 // end
902 if (v.second[i]->m_elmtVert == 0)
903 {
904 nbeg[cnt] += 1.0;
905 }
906 }
907 }
908 ++cnt;
909 }
910 Gs::Gather(nbeg, Gs::gs_add, gscomm);
911
912 // Finally determine the maximum domain id between the two vessels
913 // at an interfaces or the maximum domain id of the daughter
914 // vessels for a bifurcation or the two parent id's for a merging
915 // junction
916 Array<OneD, NekDouble> dom(VidToDomain.size(), 0.0);
917 cnt = 0;
918 for (auto &v : VidToDomain)
919 {
920 if (nvid[cnt] == 2.0)
921 {
922 // store the domain id of the two domains
923 for (size_t i = 0; i < v.second.size(); ++i)
924 {
925 dom[cnt] =
926 max(dom[cnt], (NekDouble)domId[v.second[i]->m_domain]);
927 }
928 }
929 else if (nvid[cnt] == 3.0)
930 {
931 // store the domain id of the two daughter vessels if a
932 // bifurcation or the two parent vessels if a merging
933 // junction
934
935 int val = (nbeg[cnt] == 2.0) ? 0 : 1;
936
937 for (size_t i = 0; i < v.second.size(); ++i)
938 {
939 if (v.second[i]->m_elmtVert == val)
940 {
941 dom[cnt] =
942 max(dom[cnt], (NekDouble)domId[v.second[i]->m_domain]);
943 }
944 }
945 }
946 ++cnt;
947 }
948 Gs::Gather(dom, Gs::gs_max, gscomm);
949
950 // loop over map and set up Interface information;
951 int v = 0;
952 for (auto &iter : VidToDomain)
953 {
954 ASSERTL1(nvid[v] != 1.0, "Found an interface wth only "
955 "one domain which should not happen");
956
957 if (nvid[v] == 2.0) // Vessel jump interface
958 {
959 for (size_t i = 0; i < iter.second.size(); ++i)
960 {
961 if (domId[iter.second[i]->m_domain] == dom[v])
962 {
963 iter.second[i]->m_riemannOrd = 1;
964 }
965 else
966 {
967 iter.second[i]->m_riemannOrd = 0;
968 }
969 }
970 m_vesselIntfcs.push_back(iter.second);
971 }
972 else if (nvid[v] == 3.0) // Bifurcation or Merging junction.
973 {
974 // Set up Bifurcation information
975 int val = (nbeg[v] == 2.0) ? 1 : 0;
976
977 for (size_t i = 0; i < iter.second.size(); ++i)
978 {
979 if (iter.second[i]->m_elmtVert == val)
980 {
981 iter.second[i]->m_riemannOrd = 0;
982 }
983 else
984 {
985 if (domId[iter.second[i]->m_domain] == dom[v])
986 {
987 iter.second[i]->m_riemannOrd = 2;
988 }
989 else
990 {
991 iter.second[i]->m_riemannOrd = 1;
992 }
993 }
994 }
995
996 if (nbeg[v] == 2.0)
997 {
998 m_bifurcations.push_back(iter.second);
999 }
1000 else
1001 {
1002 m_mergingJcts.push_back(iter.second);
1003 }
1004 }
1005 else
1006 {
1007 ASSERTL0(false, "Unknown junction type");
1008 }
1009 ++v; // incrementing loop over vertices
1010 }
1011
1012 // finally set up a communicator for interface data collection
1013 Array<OneD, long> tmp1(3 * nvid2dom);
1014 if (nvid2dom)
1015 {
1016 Vmath::Zero(3 * nvid2dom, tmp1, 1);
1017 }
1018
1019 cnt = 0;
1020 for (size_t n = 0; n < m_bifurcations.size(); ++n, ++cnt)
1021 {
1022 size_t vid = m_bifurcations[n][0]->m_vid;
1023 for (size_t i = 0; i < 3; ++i)
1024 {
1025 tmp1[3 * cnt + i] = (NekDouble)(3 * vid + i + 1);
1026 }
1027 }
1028
1029 for (size_t n = 0; n < m_mergingJcts.size(); ++n, ++cnt)
1030 {
1031 size_t vid = m_mergingJcts[n][0]->m_vid;
1032 for (size_t i = 0; i < 3; ++i)
1033 {
1034 tmp1[3 * cnt + i] = (NekDouble)(3 * vid + i + 1);
1035 }
1036 }
1037
1038 for (size_t n = 0; n < m_vesselIntfcs.size(); ++n, ++cnt)
1039 {
1040 size_t vid = m_vesselIntfcs[n][0]->m_vid;
1041 for (size_t i = 0; i < 3; ++i)
1042 {
1043 tmp1[3 * cnt + i] = (NekDouble)(3 * vid + i + 1);
1044 }
1045 }
1046
1047 m_intComm = Gs::Init(tmp1, m_comm->GetRowComm(), verbose);
1048}
1049
1052 map<int, LibUtilities::CommSharedPtr> &domComm)
1053{
1054
1056 Allbcs.GetBoundaryRegions();
1057
1058 /* Loop over domain and find out if we have any undefined
1059 * boundary conditions representing interfaces. If so make a
1060 * map based around vid and storing the domains that are
1061 * part of interfaces.
1062 */
1063 int numNewBc = 1;
1064 for (auto &d : m_domOrder)
1065 {
1066
1067 map<int, SpatialDomains::GeometrySharedPtr> domvids;
1068
1069 // Loop over each domain and find which vids are only touched
1070 // by one element
1071 for (auto &compIt : m_domain[d])
1072 {
1073 for (size_t j = 0; j < compIt.second->m_geomVec.size(); ++j)
1074 {
1075
1076 // get hold of vids of each segment
1077 for (size_t p = 0; p < 2; ++p)
1078 {
1080 compIt.second->m_geomVec[j]->GetVertex(p);
1081 int vid = vert->GetGlobalID();
1082
1083 // if vid has already been added remove it otherwise add
1084 // into set
1085 if (domvids.count(vid))
1086 {
1087 domvids.erase(vid);
1088 }
1089 else
1090 {
1091 domvids[vid] = vert;
1092 }
1093 }
1094 }
1095 }
1096
1097 ASSERTL1(domvids.size() == 2,
1098 "Failed to find two end points "
1099 "of a domain (domvids = " +
1100 boost::lexical_cast<std::string>(domvids.size()) + ")");
1101
1102 size_t nprocs = domComm[d]->GetSize();
1103
1104 if (nprocs > 1) // Remove parallel interfaces
1105 {
1106 size_t rank = domComm[d]->GetRank();
1107 Array<OneD, int> nvids(nprocs, 0);
1108 nvids[rank] = domvids.size();
1109 domComm[d]->AllReduce(nvids, LibUtilities::ReduceSum);
1110
1111 size_t totvids = Vmath::Vsum(nprocs, nvids, 1);
1112
1113 Array<OneD, int> locids(totvids, -1);
1114
1115 size_t cnt = 0;
1116 for (size_t i = 0; i < rank; ++i)
1117 {
1118 cnt += nvids[i];
1119 }
1120
1121 for (auto &i : domvids)
1122 {
1123 locids[cnt++] = i.first;
1124 }
1125
1126 // collect lcoal vids on this domain on all processor
1127 domComm[d]->AllReduce(locids, LibUtilities::ReduceMax);
1128
1129 set<int> chkvids;
1130 for (size_t i = 0; i < totvids; ++i)
1131 {
1132 if (chkvids.count(locids[i]))
1133 {
1134 // if this id is on local domain then remove
1135 if (domvids.count(locids[i]))
1136 {
1137 domvids.erase(locids[i]);
1138 }
1139 }
1140 else
1141 {
1142 chkvids.insert(locids[i]);
1143 }
1144 }
1145 }
1146
1147 // Finally check if remaining domvids are already declared, if
1148 // not add them
1149
1150 // Erase from domvids any existing BCs
1151 for (auto &it : bregions)
1152 {
1153 for (auto &bregionIt : *it.second)
1154 {
1155 // can assume that all regions only contain one point in 1D
1156 // Really do not need loop above
1157 int id = bregionIt.second->m_geomVec[0]->GetGlobalID();
1158
1159 if (domvids.count(id))
1160 {
1161 domvids.erase(id);
1162 }
1163 }
1164 }
1165
1166 // Finally add interface conditions to Allbcs
1167 std::vector<std::string> variables = m_session->GetVariables();
1168
1169 // determine the maxmimum bregion index
1170 int maxbregind = 0;
1171 for (auto &b : bregions)
1172 {
1173 maxbregind = (maxbregind > b.first) ? maxbregind : b.first;
1174 }
1175
1176 for (auto &b : domvids)
1177 {
1180 SpatialDomains::BoundaryRegion>::AllocateSharedPtr());
1181
1182 // Set up Composite (GemetryVector) to contain vertex and put into
1183 // bRegion
1186 gvec->m_geomVec.push_back(b.second);
1187 (*breg)[b.first] = gvec;
1188
1189 Allbcs.AddBoundaryRegions(maxbregind + numNewBc, breg);
1190
1193 SpatialDomains::BoundaryConditionMap>::AllocateSharedPtr();
1194
1195 std::string userDefined = "Interface";
1196 // Set up just boundary condition for this variable.
1197 SpatialDomains::BoundaryConditionShPtr DirichletInterface(
1199 AllocateSharedPtr(m_session, "0", userDefined));
1200
1201 for (size_t i = 0; i < variables.size(); ++i)
1202 {
1203 (*bCondition)[variables[i]] = DirichletInterface;
1204 }
1205
1206 Allbcs.AddBoundaryConditions(maxbregind + numNewBc, bCondition);
1207 ++numNewBc;
1208 }
1209 }
1210}
1211
1212/**
1213 * Initialisation routine for multiple subdomain case. Sets the
1214 * initial conditions for all arterial subdomains read from the
1215 * inputfile. Sets the material properties and the A_0 area for
1216 * all subdomains and fills the domain-linking boundary
1217 * conditions with the initial values of their domain.
1218 */
1219void PulseWaveSystem::v_DoInitialise(bool dumpInitialConditions)
1220{
1221 boost::ignore_unused(dumpInitialConditions);
1222
1223 if (m_session->GetComm()->GetRank() == 0)
1224 {
1225 cout << "Initial Conditions: " << endl;
1226 }
1227
1228 /* Loop over all subdomains to initialize all with the Initial
1229 * Conditions read from the inputfile*/
1230 int omega = 0;
1231 for (auto &d : m_domOrder)
1232 {
1233 for (size_t i = 0; i < 2; ++i)
1234 {
1235 m_fields[i] = m_vessels[m_nVariables * omega + i];
1236 }
1237
1238 if (m_session->GetComm()->GetRank() == 0)
1239 {
1240 cout << "Subdomain: " << omega << endl;
1241 }
1242
1243 SetInitialConditions(0.0, false, d);
1244 omega++;
1245 }
1246
1247 // Reset 2 variables to first vessels
1248 for (size_t i = 0; i < 2; ++i)
1249 {
1250 m_fields[i] = m_vessels[i];
1251 }
1252}
1253
1254/**
1255 * NEEDS Updating:
1256 *
1257 * DoSolve routine for PulseWavePropagation with multiple
1258 * subdomains taken from UnsteadySystem and modified for
1259 * multidomain case. Initialises the time integration scheme (as
1260 * specified in the session file), and perform the time
1261 * integration. Within the timestepping loop the following is
1262 * done: 1. Link all arterial segments according to the network
1263 * structure, solve the Riemann problem between different
1264 * arterial segments and assign the values to the boundary
1265 * conditions (LinkSubdomains) 2. Every arterial segment is
1266 * solved independently for this timestep. This is done by handing
1267 * the solution vector \f$ \mathbf{u} \f$ and the right hand side
1268 * m_ode, which is the PulseWavePropagation class in this example
1269 * over to the time integration scheme
1270 */
1272{
1273 // NekDouble IntegrationTime = 0.0;
1274 size_t i;
1275 int n;
1276 int nchk = 1;
1277
1279
1280 for (i = 0; i < m_nVariables; ++i)
1281 {
1282 fields[i] = m_vessels[i]->UpdatePhys();
1283 m_fields[i]->SetPhysState(false);
1284 }
1285
1286 m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
1287
1288 // Time loop
1289 for (n = 0; n < m_steps; ++n)
1290 {
1291 LibUtilities::Timer time_v_IntStep;
1292 time_v_IntStep.Start();
1293 fields = m_intScheme->TimeIntegrate(n, m_timestep);
1294 m_time += m_timestep;
1295 time_v_IntStep.Stop();
1296 time_v_IntStep.AccumulateRegion("PulseWaveSystem::TimeIntegrationStep",
1297 0);
1298 // IntegrationTime += timer.TimePerTest(1);
1299
1300 // Write out status information.
1301 if (m_infosteps && !((n + 1) % m_infosteps) &&
1302 m_session->GetComm()->GetRank() == 0)
1303 {
1304 cout << "Steps: " << n + 1 << "\t Time: " << m_time
1305 << "\t Time-step: " << m_timestep << "\t" << endl;
1306 }
1307
1308 // Transform data if needed
1309 if (m_checksteps && !((n + 1) % m_checksteps))
1310 {
1311 for (i = 0; i < m_nVariables; ++i)
1312 {
1313 size_t cnt = 0;
1314 for (size_t omega = 0; omega < m_nDomains; omega++)
1315 {
1316 m_vessels[omega * m_nVariables + i]->FwdTrans(
1317 fields[i] + cnt,
1318 m_vessels[omega * m_nVariables + i]->UpdateCoeffs());
1319 cnt += m_vessels[omega * m_nVariables + i]->GetTotPoints();
1320 }
1321 }
1322 CheckPoint_Output(nchk++);
1323 }
1324
1325 } // end of timeintegration
1326
1327 // Copy Array To Vessel Phys Fields
1328 for (size_t i = 0; i < m_nVariables; ++i)
1329 {
1330 Vmath::Vcopy(fields[i].size(), fields[i], 1, m_vessels[i]->UpdatePhys(),
1331 1);
1332 }
1333
1334 /** if(m_session->GetComm()->GetRank() == 0)
1335 {
1336 cout << "Time-integration timing: " << IntegrationTime << " s" <<
1337 endl
1338 << endl;
1339 } **/
1340}
1341
1344 const Array<OneD, const Array<OneD, NekDouble>> &fields, NekDouble &A,
1345 NekDouble &u, NekDouble &beta, NekDouble &A_0, NekDouble &alpha)
1346{
1347 LibUtilities::Timer time_FillDataFromInterfacePoint;
1348 time_FillDataFromInterfacePoint.Start();
1349
1350 int omega = I->m_domain;
1351 int traceId = I->m_traceId;
1352 int eid = I->m_elmt;
1353 int vert = I->m_elmtVert;
1354 int vesselID = omega * m_nVariables;
1355 int phys_offset = m_vessels[vesselID]->GetPhys_Offset(eid);
1358
1359 m_vessels[vesselID]->GetExp(eid)->GetTracePhysVals(
1360 vert, dumExp, fields[0] + m_fieldPhysOffset[omega] + phys_offset, Tmp);
1361 A = Tmp[0];
1362 m_vessels[vesselID]->GetExp(eid)->GetTracePhysVals(
1363 vert, dumExp, fields[1] + m_fieldPhysOffset[omega] + phys_offset, Tmp);
1364 u = Tmp[0];
1365
1366 beta = m_beta_trace[omega][traceId];
1367 A_0 = m_A_0_trace[omega][traceId];
1368 alpha = m_alpha_trace[omega][traceId];
1369 time_FillDataFromInterfacePoint.Stop();
1370 time_FillDataFromInterfacePoint.AccumulateRegion(
1371 "PulseWaveSystem::FillDataFromInterfacePoint", 2);
1372}
1373
1375 const Array<OneD, const Array<OneD, NekDouble>> &fields)
1376{
1377 LibUtilities::Timer time_EnforceInterfaceConditions;
1378 time_EnforceInterfaceConditions.Start();
1379 int dom, bcpos;
1380
1381 int totif =
1382 m_bifurcations.size() + m_mergingJcts.size() + m_vesselIntfcs.size();
1383
1384 Array<OneD, NekDouble> Aut, Au(3 * totif, 0.0);
1385 Array<OneD, NekDouble> uut, uu(3 * totif, 0.0);
1386 Array<OneD, NekDouble> betat, beta(3 * totif, 0.0);
1387 Array<OneD, NekDouble> A_0t, A_0(3 * totif, 0.0);
1388 Array<OneD, NekDouble> alphat, alpha(3 * totif, 0.0);
1389
1390 // Bifurcations Data:
1391 size_t cnt = 0;
1392 for (size_t n = 0; n < m_bifurcations.size(); ++n, ++cnt)
1393 {
1394 for (size_t i = 0; i < m_bifurcations[n].size(); ++i)
1395 {
1396 size_t l = m_bifurcations[n][i]->m_riemannOrd;
1398 m_bifurcations[n][i], fields, Au[l + 3 * cnt], uu[l + 3 * cnt],
1399 beta[l + 3 * cnt], A_0[l + 3 * cnt], alpha[l + 3 * cnt]);
1400 }
1401 }
1402
1403 // Enforce Merging vessles Data:
1404 for (size_t n = 0; n < m_mergingJcts.size(); ++n, ++cnt)
1405 {
1406 // Merged vessel
1407 for (size_t i = 0; i < m_mergingJcts.size(); ++i)
1408 {
1409 int l = m_mergingJcts[n][i]->m_riemannOrd;
1411 m_mergingJcts[n][i], fields, Au[l + 3 * cnt], uu[l + 3 * cnt],
1412 beta[l + 3 * cnt], A_0[l + 3 * cnt], alpha[l + 3 * cnt]);
1413 }
1414 }
1415
1416 // Enforce interface junction between two vessesls
1417 for (size_t n = 0; n < m_vesselIntfcs.size(); ++n, ++cnt)
1418 {
1419 for (size_t i = 0; i < m_vesselIntfcs[n].size(); ++i)
1420 {
1421 int l = m_vesselIntfcs[n][i]->m_riemannOrd;
1423 m_vesselIntfcs[n][i], fields, Au[l + 3 * cnt], uu[l + 3 * cnt],
1424 beta[l + 3 * cnt], A_0[l + 3 * cnt], alpha[l + 3 * cnt]);
1425 }
1426 }
1427
1428 // Gather data if running in parallel
1434
1435 // Enforce Bifurcations:
1436 cnt = 0;
1437 for (size_t n = 0; n < m_bifurcations.size(); ++n, ++cnt)
1438 {
1439 // Solve the Riemann problem for a bifurcation
1440 BifurcationRiemann(Aut = Au + 3 * cnt, uut = uu + 3 * cnt,
1441 betat = beta + 3 * cnt, A_0t = A_0 + 3 * cnt,
1442 alphat = alpha + 3 * cnt);
1443
1444 // Store the values into the right positions:
1445 for (size_t i = 0; i < m_bifurcations[n].size(); ++i)
1446 {
1447 dom = m_bifurcations[n][i]->m_domain;
1448 bcpos = m_bifurcations[n][i]->m_bcPosition;
1449 int l = m_bifurcations[n][i]->m_riemannOrd;
1450 m_vessels[dom * m_nVariables]
1451 ->UpdateBndCondExpansion(bcpos)
1452 ->UpdatePhys()[0] = Au[l + 3 * cnt];
1453 m_vessels[dom * m_nVariables + 1]
1454 ->UpdateBndCondExpansion(bcpos)
1455 ->UpdatePhys()[0] = uu[l + 3 * cnt];
1456 }
1457 }
1458
1459 // Enforce Merging vessles;
1460 for (size_t n = 0; n < m_mergingJcts.size(); ++n, ++cnt)
1461 {
1462 // Solve the Riemann problem for a merging vessel
1463 MergingRiemann(Aut = Au + 3 * cnt, uut = uu + 3 * cnt,
1464 betat = beta + 3 * cnt, A_0t = A_0 + 3 * cnt,
1465 alphat = alpha + 3 * cnt);
1466
1467 // Store the values into the right positions:
1468 for (size_t i = 0; i < m_mergingJcts[n].size(); ++i)
1469 {
1470 int dom = m_mergingJcts[n][i]->m_domain;
1471 int bcpos = m_mergingJcts[n][i]->m_bcPosition;
1472 int l = m_mergingJcts[n][i]->m_riemannOrd;
1473 m_vessels[dom * m_nVariables]
1474 ->UpdateBndCondExpansion(bcpos)
1475 ->UpdatePhys()[0] = Au[l + 3 * cnt];
1476 m_vessels[dom * m_nVariables + 1]
1477 ->UpdateBndCondExpansion(bcpos)
1478 ->UpdatePhys()[0] = uu[l + 3 * cnt];
1479 }
1480 }
1481
1482 // Enforce interface junction between two vessesls
1483 for (size_t n = 0; n < m_vesselIntfcs.size(); ++n, ++cnt)
1484 {
1485 InterfaceRiemann(Aut = Au + 3 * cnt, uut = uu + 3 * cnt,
1486 betat = beta + 3 * cnt, A_0t = A_0 + 3 * cnt,
1487 alphat = alpha + 3 * cnt);
1488
1489 // Store the values into the right positions:
1490 for (size_t i = 0; i < m_vesselIntfcs[n].size(); ++i)
1491 {
1492 int dom = m_vesselIntfcs[n][i]->m_domain;
1493 int bcpos = m_vesselIntfcs[n][i]->m_bcPosition;
1494 int l = m_vesselIntfcs[n][i]->m_riemannOrd;
1495 m_vessels[dom * m_nVariables]
1496 ->UpdateBndCondExpansion(bcpos)
1497 ->UpdatePhys()[0] = Au[l + 3 * cnt];
1498 m_vessels[dom * m_nVariables + 1]
1499 ->UpdateBndCondExpansion(bcpos)
1500 ->UpdatePhys()[0] = uu[l + 3 * cnt];
1501 }
1502 }
1503 time_EnforceInterfaceConditions.Stop();
1504 time_EnforceInterfaceConditions.AccumulateRegion(
1505 "PulseWaveSystem::EnforceInterfaceConditions", 1);
1506}
1507
1508/**
1509 * Solves the Riemann problem at a bifurcation by assuming
1510 * subsonic flow at both sides of the boundary and by
1511 * applying conservation of mass and continuity of the total
1512 * pressure \f$ \frac{p}{rho} + \frac{u^{2}}{2}. \f$ The
1513 * other 3 missing equations come from the characteristic
1514 * variables. For further information see "Pulse
1515 * WavePropagation in the human vascular system" Section
1516 * 3.4.4
1517 */
1523{
1524
1525 NekDouble rho = m_rho;
1527 Array<OneD, NekDouble> P_Au(3);
1528 Array<OneD, NekDouble> W_Au(3);
1529 NekMatrix<NekDouble> invJ(6, 6);
1532
1533 int proceed = 1;
1534 int iter = 0;
1535 int MAX_ITER = 100;
1536
1537 // Forward characteristic
1538 m_pressureArea->GetW1(W[0], uu[0], beta[0], Au[0], A_0[0], alpha[0]);
1539
1540 // Backward characteristics
1541 for (size_t i = 1; i < 3; ++i)
1542 {
1543 m_pressureArea->GetW2(W[i], uu[i], beta[i], Au[i], A_0[i], alpha[i]);
1544 }
1545
1546 // Tolerances for the algorithm
1547 NekDouble Tol = 1.0E-10;
1548
1549 // Newton Iteration
1550 while ((proceed) && (iter < MAX_ITER))
1551 {
1552 LibUtilities::Timer time_BifurcationRiemann;
1553 time_BifurcationRiemann.Start();
1554 iter += 1;
1555
1556 /*
1557 * We solve the six constraint equations via a multivariate Newton
1558 * iteration. Equations are:
1559 * 1. Forward characteristic: W1(A_L, U_L) = W1(Au_L, Uu_L)
1560 * 2. Backward characteristic 1: W2(A_R1, U_R1) = W2(Au_R1, Uu_R1)
1561 * 3. Backward characteristic 2: W2(A_R2, U_R2) = W2(Au_R2, Uu_R2)
1562 * 4. Conservation of mass: Au_L * Uu_L = Au_R1 * Uu_R1 +
1563 * Au_R2 * Uu_R2
1564 * 5. Continuity of total pressure 1: rho * Uu_L * Uu_L / 2 + p(Au_L)
1565 * = rho * Uu_R1 * Uu_R1 / 2 + p(Au_R1)
1566 * 6. Continuity of total pressure 2: rho * Uu_L * Uu_L / 2 + p(Au_L)
1567 * = rho * Uu_R2 * Uu_R2 / 2 + p(Au_R2)
1568 */
1569 for (size_t i = 0; i < 3; ++i)
1570 {
1571 m_pressureArea->GetPressure(P_Au[i], beta[i], Au[i], A_0[i], 0, 0,
1572 alpha[i]);
1573 }
1574
1575 m_pressureArea->GetW1(W_Au[0], uu[0], beta[0], Au[0], A_0[0], alpha[0]);
1576 m_pressureArea->GetW2(W_Au[1], uu[1], beta[1], Au[1], A_0[1], alpha[1]);
1577 m_pressureArea->GetW2(W_Au[2], uu[2], beta[2], Au[2], A_0[2], alpha[2]);
1578
1579 // Constraint equations set to zero
1580 f[0] = W_Au[0] - W[0];
1581 f[1] = W_Au[1] - W[1];
1582 f[2] = W_Au[2] - W[2];
1583 f[3] = Au[0] * uu[0] - Au[1] * uu[1] - Au[2] * uu[2];
1584 f[4] = uu[0] * uu[0] + 2 * P_Au[0] / rho - uu[1] * uu[1] -
1585 2 * P_Au[1] / rho;
1586 f[5] = uu[0] * uu[0] + 2 * P_Au[0] / rho - uu[2] * uu[2] -
1587 2 * P_Au[2] / rho;
1588
1589 // Inverse Jacobian for x + dx = x - J^(-1) * f(x)
1590 m_pressureArea->GetJacobianInverse(invJ, Au, uu, beta, A_0, alpha,
1591 "Bifurcation");
1592
1593 Multiply(dx, invJ, f);
1594
1595 // Update the solution: x_new = x_old - dx
1596 for (int i = 0; i < 3; ++i)
1597 {
1598 uu[i] -= dx[i];
1599 Au[i] -= dx[i + 3];
1600 }
1601
1602 // Check if the error of the solution is smaller than Tol
1603 if (Dot(dx, dx) < Tol)
1604 {
1605 proceed = 0;
1606 }
1607
1608 // Check if solver converges
1609 if (iter >= MAX_ITER)
1610 {
1611 ASSERTL0(false, "Riemann solver for Bifurcation did not converge.");
1612 }
1613 time_BifurcationRiemann.Stop();
1614 time_BifurcationRiemann.AccumulateRegion(
1615 "PulseWaveSystem::Bifurcation Riemann", 2);
1616 }
1617}
1618
1619/**
1620 * Solves the Riemann problem at an merging flow condition by
1621 * assuming subsonic flow at both sides of the boundary and by
1622 * applying conservation of mass and continuity of the total
1623 * pressure \f$ \frac{p}{rho} + \frac{u^{2}}{2}. \f$ The other 3
1624 * missing equations come from the characteristic variables. For
1625 * further information see "Pulse WavePropagation in the human
1626 * vascular system" Section 3.4.4
1627 */
1633{
1634
1635 NekDouble rho = m_rho;
1637 Array<OneD, NekDouble> W_Au(3);
1638 Array<OneD, NekDouble> P_Au(3);
1639 NekMatrix<NekDouble> invJ(6, 6);
1642
1643 int proceed = 1;
1644 int iter = 0;
1645 int MAX_ITER = 15;
1646
1647 // Backward characteristic
1648 m_pressureArea->GetW2(W[0], uu[0], beta[0], Au[0], A_0[0], alpha[0]);
1649
1650 // Forward characteristics
1651 for (size_t i = 1; i < 3; ++i)
1652 {
1653 m_pressureArea->GetW1(W[i], uu[i], beta[i], Au[i], A_0[i], alpha[i]);
1654 }
1655
1656 // Tolerances for the algorithm
1657 NekDouble Tol = 1.0E-10;
1658
1659 // Newton Iteration
1660 while ((proceed) && (iter < MAX_ITER))
1661 {
1662 LibUtilities::Timer time_MergingRiemann;
1663 time_MergingRiemann.Start();
1664 iter += 1;
1665
1666 /*
1667 * We solve the six constraint equations via a multivariate Newton
1668 * iteration. Equations are:
1669 * 1. Backward characteristic: W2(A_R, U_R) = W1(Au_R, Uu_R)
1670 * 2. Forward characteristic 1: W1(A_L1, U_L1) = W1(Au_L1,
1671 * Uu_R1)
1672 * 3. Forward characteristic 2: W1(A_L2, U_L2) = W1(Au_L2,
1673 * Uu_L2)
1674 * 4. Conservation of mass: Au_R * Uu_R = Au_L1 * Uu_L1 +
1675 * Au_L2 * Uu_L2
1676 * 5. Continuity of total pressure 1: rho * Uu_R * Uu_R / 2 + p(Au_R)
1677 * = rho * Uu_L1 * Uu_L1 / 2 + p(Au_L1)
1678 * 6. Continuity of total pressure 2: rho * Uu_R * Uu_R / 2 + p(Au_R)
1679 * = rho * Uu_L2 * Uu_L2 / 2 + p(Au_L2)
1680 */
1681 for (size_t i = 0; i < 3; ++i)
1682 {
1683 m_pressureArea->GetPressure(P_Au[i], beta[i], Au[i], A_0[i], 0, 0,
1684 alpha[i]);
1685 }
1686
1687 m_pressureArea->GetW2(W_Au[0], uu[0], beta[0], Au[0], A_0[0], alpha[0]);
1688 m_pressureArea->GetW1(W_Au[1], uu[1], beta[1], Au[1], A_0[1], alpha[1]);
1689 m_pressureArea->GetW1(W_Au[2], uu[2], beta[2], Au[2], A_0[2], alpha[2]);
1690
1691 // Constraint equations set to zero
1692 f[0] = W_Au[0] - W[0];
1693 f[1] = W_Au[1] - W[1];
1694 f[2] = W_Au[2] - W[2];
1695 f[3] = Au[0] * uu[0] - Au[1] * uu[1] - Au[2] * uu[2];
1696 f[4] = uu[0] * uu[0] + 2 * P_Au[0] / rho - uu[1] * uu[1] -
1697 2 * P_Au[1] / rho;
1698 f[5] = uu[0] * uu[0] + 2 * P_Au[0] / rho - uu[2] * uu[2] -
1699 2 * P_Au[2] / rho;
1700
1701 // Inverse Jacobian for x + dx = x - J^(-1) * f(x)
1702 m_pressureArea->GetJacobianInverse(invJ, Au, uu, beta, A_0, alpha,
1703 "Merge");
1704
1705 Multiply(dx, invJ, f);
1706
1707 // Update the solution: x_new = x_old - dx
1708 for (size_t i = 0; i < 3; ++i)
1709 {
1710 uu[i] -= dx[i];
1711 Au[i] -= dx[i + 3];
1712 }
1713
1714 // Check if the error of the solution is smaller than Tol
1715 if (Dot(dx, dx) < Tol)
1716 {
1717 proceed = 0;
1718 }
1719
1720 // Check if solver converges
1721 if (iter >= MAX_ITER)
1722 {
1723 ASSERTL0(false, "Riemann solver for Merging Flow did not converge");
1724 }
1725 time_MergingRiemann.Stop();
1726 time_MergingRiemann.AccumulateRegion("PulseWaveSystem::MergingRiemann",
1727 2);
1728 }
1729}
1730
1731/**
1732 * Solves the Riemann problem at an interdomain
1733 * junction/Interface by assuming subsonic flow at both sides of
1734 * the boundary and by applying conservation of mass and
1735 * continuity of the total pressure \f$ \frac{p}{rho} +
1736 * \frac{u^{2}}{2}. \f$ The other 2 missing equations come from
1737 * the characteristic variables. For further information see
1738 * "Pulse WavePropagation in the human vascular system" Section
1739 * 3.4.
1740 */
1746{
1747
1748 NekDouble rho = m_rho;
1750 Array<OneD, NekDouble> W_Au(2);
1751 Array<OneD, NekDouble> P_Au(2);
1752 NekMatrix<NekDouble> invJ(4, 4);
1755
1756 int proceed = 1;
1757 int iter = 0;
1758 int MAX_ITER = 15;
1759 NekDouble Tol = 1.0E-10;
1760
1761 // Forward and backward characteristics
1762 m_pressureArea->GetW1(W[0], uu[0], beta[0], Au[0], A_0[0], alpha[0]);
1763 m_pressureArea->GetW2(W[1], uu[1], beta[1], Au[1], A_0[1], alpha[1]);
1764
1765 while ((proceed) && (iter < MAX_ITER))
1766 {
1767 LibUtilities::Timer time_InterfaceRiemann;
1768 time_InterfaceRiemann.Start();
1769 iter += 1;
1770
1771 /*
1772 * We solve the four constraint equations via a multivariate Newton
1773 * iteration. Equations are:
1774 * 1. Forward characteristic: W1(A_L, U_L) = W1(Au_L, Uu_L)
1775 * 2. Backward characteristic: W2(A_R, U_R) = W2(Au_R, Uu_R)
1776 * 3. Conservation of mass: Au_L * Uu_L = Au_R * Uu_R
1777 * 4. Continuity of total pressure: rho * Uu_L * Uu_L / 2 + p(Au_L) =
1778 * rho * Uu_R * Uu_R / 2 + p(Au_R)
1779 */
1780 for (size_t i = 0; i < 2; ++i)
1781 {
1782 m_pressureArea->GetPressure(P_Au[i], beta[i], Au[i], A_0[i], 0, 0,
1783 alpha[i]);
1784 }
1785
1786 m_pressureArea->GetW1(W_Au[0], uu[0], beta[0], Au[0], A_0[0], alpha[0]);
1787 m_pressureArea->GetW2(W_Au[1], uu[1], beta[1], Au[1], A_0[1], alpha[1]);
1788
1789 // Constraint equations set to zero
1790 f[0] = W_Au[0] - W[0];
1791 f[1] = W_Au[1] - W[1];
1792 f[2] = Au[0] * uu[0] - Au[1] * uu[1];
1793 f[3] = uu[0] * uu[0] + 2 * P_Au[0] / rho - uu[1] * uu[1] -
1794 2 * P_Au[1] / rho;
1795
1796 // Inverse Jacobian for x + dx = x - J^(-1) * f(x)
1797 m_pressureArea->GetJacobianInverse(invJ, Au, uu, beta, A_0, alpha,
1798 "Interface");
1799
1800 Multiply(dx, invJ, f);
1801
1802 // Update solution: x_new = x_old - dx
1803 for (size_t i = 0; i < 2; ++i)
1804 {
1805 uu[i] -= dx[i];
1806 Au[i] -= dx[i + 2];
1807 }
1808
1809 // Check if the error of the solution is smaller than Tol.
1810 if (Dot(dx, dx) < Tol)
1811 {
1812 proceed = 0;
1813 }
1814 time_InterfaceRiemann.Stop();
1815 time_InterfaceRiemann.AccumulateRegion(
1816 "PulseWaveSystem::InterfaceRiemann", 2);
1817 }
1818
1819 if (iter >= MAX_ITER)
1820 {
1821 ASSERTL0(false, "Riemann solver for Interface did not converge");
1822 }
1823}
1824
1825/**
1826 * Writes the .fld file at the end of the simulation. Similar to the normal
1827 * v_Output however the Multidomain output has to be prepared.
1828 */
1830{
1831 /**
1832 * Write the field data to file. The file is named according to the session
1833 * name with the extension .fld appended.
1834 */
1835 std::string outname = m_sessionName + ".fld";
1836
1837 WriteVessels(outname);
1838}
1839
1840/**
1841 * Writes the .fld file at the end of the simulation. Similar to the normal
1842 * v_Output however the Multidomain output has to be prepared.
1843 */
1845{
1846 std::stringstream outname;
1847 outname << m_sessionName << "_" << n << ".chk";
1848
1849 WriteVessels(outname.str());
1850}
1851
1852/**
1853 * Writes the field data to a file with the given filename.
1854 * @param outname Filename to write to.
1855 */
1856void PulseWaveSystem::WriteVessels(const std::string &outname)
1857{
1858 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1859 std::vector<std::string> variables = m_session->GetVariables();
1860
1861 for (size_t n = 0; n < m_nDomains; ++n)
1862 {
1863 m_vessels[n * m_nVariables]->GetFieldDefinitions(FieldDef);
1864 }
1865
1866 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
1867
1868 size_t nFieldDefPerDomain = FieldDef.size() / m_nDomains;
1869 size_t cnt;
1870
1871 // Copy Data into FieldData and set variable
1872 for (size_t n = 0; n < m_nDomains; ++n)
1873 {
1874 // Outputs area and velocity
1875 for (size_t j = 0; j < m_nVariables; ++j)
1876 {
1877 for (size_t i = 0; i < nFieldDefPerDomain; ++i)
1878 {
1879 cnt = n * nFieldDefPerDomain + i;
1880
1881 FieldDef[cnt]->m_fields.push_back(variables[j]);
1882
1883 m_vessels[n * m_nVariables]->AppendFieldData(
1884 FieldDef[cnt], FieldData[cnt],
1885 m_vessels[n * m_nVariables + j]->UpdateCoeffs());
1886 }
1887 }
1888
1889 // Outputs pressure
1891
1892 m_vessels[n * m_nVariables]->FwdTransLocalElmt(m_pressure[n], PFwd);
1893
1894 for (size_t i = 0; i < nFieldDefPerDomain; ++i)
1895 {
1896 cnt = n * nFieldDefPerDomain + i;
1897
1898 FieldDef[cnt]->m_fields.push_back("P");
1899
1900 m_vessels[n * m_nVariables]->AppendFieldData(FieldDef[cnt],
1901 FieldData[cnt], PFwd);
1902 }
1903
1904 if (extraFields)
1905 {
1912
1913 m_vessels[n * m_nVariables]->FwdTransLocalElmt(m_PWV[n], PWVFwd);
1914 m_vessels[n * m_nVariables]->FwdTransLocalElmt(m_W1[n], W1Fwd);
1915 m_vessels[n * m_nVariables]->FwdTransLocalElmt(m_W2[n], W2Fwd);
1916
1917 for (size_t i = 0; i < nFieldDefPerDomain; ++i)
1918 {
1919 cnt = n * nFieldDefPerDomain + i;
1920
1921 FieldDef[cnt]->m_fields.push_back("c");
1922 FieldDef[cnt]->m_fields.push_back("W1");
1923 FieldDef[cnt]->m_fields.push_back("W2");
1924
1925 m_vessels[n * m_nVariables]->AppendFieldData(
1926 FieldDef[cnt], FieldData[cnt], PWVFwd);
1927 m_vessels[n * m_nVariables]->AppendFieldData(
1928 FieldDef[cnt], FieldData[cnt], W1Fwd);
1929 m_vessels[n * m_nVariables]->AppendFieldData(
1930 FieldDef[cnt], FieldData[cnt], W2Fwd);
1931 }
1932 }
1933 }
1934
1935 // Update time in field info if required
1936 if (m_fieldMetaDataMap.find("Time") != m_fieldMetaDataMap.end())
1937 {
1938 m_fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(m_time);
1939 }
1940
1941 m_fld->Write(outname, FieldDef, FieldData, m_fieldMetaDataMap);
1942}
1943
1944/* Compute the error in the L2-norm
1945 * @param field The field to compare.
1946 * @param exactsoln The exact solution to compare with.
1947 * @param Normalised Normalise L2-error.
1948 * @returns Error in the L2-norm.
1949 */
1951 const Array<OneD, NekDouble> &exact,
1952 bool Normalised)
1953{
1954 boost::ignore_unused(exact);
1955
1956 NekDouble L2error = 0.0;
1957 NekDouble L2error_dom;
1958 NekDouble Vol = 0.0;
1959
1960 if (m_NumQuadPointsError == 0)
1961 {
1962 for (size_t omega = 0; omega < m_nDomains; omega++)
1963 {
1964 size_t vesselid = field + omega * m_nVariables;
1965
1966 // since exactsoln is passed for just the first field size we need
1967 // to reset it each domain
1968 Array<OneD, NekDouble> exactsoln(m_vessels[vesselid]->GetNpoints());
1969 m_fields[field] = m_vessels[omega * m_nVariables];
1970 EvaluateExactSolution(field, exactsoln, m_time);
1971
1972 if (m_vessels[vesselid]->GetPhysState() == false)
1973 {
1974 m_vessels[vesselid]->BwdTrans(
1975 m_vessels[vesselid]->GetCoeffs(),
1976 m_vessels[vesselid]->UpdatePhys());
1977 }
1978
1979 if (exactsoln.size())
1980 {
1981 L2error_dom = m_vessels[vesselid]->L2(
1982 m_vessels[vesselid]->GetPhys(), exactsoln);
1983 }
1984 else if (m_session->DefinesFunction("ExactSolution"))
1985 {
1986 Array<OneD, NekDouble> exactsoln(
1987 m_vessels[vesselid]->GetNpoints());
1988
1990 m_session->GetFunction("ExactSolution", field, omega);
1991 GetFunction("ExactSolution")
1992 ->Evaluate(m_session->GetVariable(field), exactsoln,
1993 m_time);
1994
1995 L2error_dom = m_vessels[vesselid]->L2(
1996 m_vessels[vesselid]->GetPhys(), exactsoln);
1997 }
1998 else
1999 {
2000 L2error_dom =
2001 m_vessels[vesselid]->L2(m_vessels[vesselid]->GetPhys());
2002 }
2003
2004 if (m_vessels[vesselid]->GetComm()->GetRank())
2005 {
2006 // ensure domains are only summed on local root of
2007 // domain
2008 L2error_dom = 0.0;
2009 }
2010
2011 L2error += L2error_dom * L2error_dom;
2012
2013 if (Normalised == true)
2014 {
2016 1.0);
2017
2018 Vol += m_vessels[vesselid]->Integral(one);
2019 }
2020 }
2021 }
2022 else
2023 {
2024 ASSERTL0(false, "Not set up");
2025 }
2026
2027 m_comm->AllReduce(L2error, LibUtilities::ReduceSum);
2028
2029 if (Normalised == true)
2030 {
2031 m_comm->AllReduce(Vol, LibUtilities::ReduceSum);
2032
2033 L2error = sqrt(L2error / Vol);
2034 }
2035 else
2036 {
2037 L2error = sqrt(L2error);
2038 }
2039
2040 return L2error;
2041}
2042
2043/**
2044 * Compute the error in the L_inf-norm
2045 * @param field The field to compare.
2046 * @param exactsoln The exact solution to compare with.
2047 * @returns Error in the L_inft-norm.
2048 */
2050 const Array<OneD, NekDouble> &exact)
2051{
2052 boost::ignore_unused(exact);
2053
2054 NekDouble LinferrorDom, Linferror = -1.0;
2055
2056 for (size_t omega = 0; omega < m_nDomains; ++omega)
2057 {
2058 size_t vesselid = field + omega * m_nVariables;
2059
2060 // since exactsoln is passed for just the first field size we need
2061 // to reset it each domain
2062 Array<OneD, NekDouble> exactsoln(m_vessels[vesselid]->GetNpoints());
2063 m_fields[field] = m_vessels[omega * m_nVariables];
2064 EvaluateExactSolution(field, exactsoln, m_time);
2065
2066 if (m_NumQuadPointsError == 0)
2067 {
2068 if (m_vessels[vesselid]->GetPhysState() == false)
2069 {
2070 m_vessels[vesselid]->BwdTrans(
2071 m_vessels[vesselid]->GetCoeffs(),
2072 m_vessels[vesselid]->UpdatePhys());
2073 }
2074
2075 if (exactsoln.size())
2076 {
2077 LinferrorDom = m_vessels[vesselid]->Linf(
2078 m_vessels[vesselid]->GetPhys(), exactsoln);
2079 }
2080 else if (m_session->DefinesFunction("ExactSolution"))
2081 {
2082 Array<OneD, NekDouble> exactsoln(
2083 m_vessels[vesselid]->GetNpoints());
2084
2085 GetFunction("ExactSolution")
2086 ->Evaluate(m_session->GetVariable(field), exactsoln,
2087 m_time);
2088
2089 LinferrorDom = m_vessels[vesselid]->Linf(
2090 m_vessels[vesselid]->GetPhys(), exactsoln);
2091 }
2092 else
2093 {
2094 LinferrorDom = 0.0;
2095 }
2096
2097 Linferror = (Linferror > LinferrorDom) ? Linferror : LinferrorDom;
2098 }
2099 else
2100 {
2101 ASSERTL0(false, "ErrorExtraPoints not allowed for this solver");
2102 }
2103 }
2104 m_comm->AllReduce(Linferror, LibUtilities::ReduceMax);
2105 return Linferror;
2106}
2107
2108} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:72
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, Array< OneD, NekDouble > > m_A_0
Array< OneD, int > m_fieldPhysOffset
virtual void v_DoInitialise(bool dumpInitialConditions=false) override
Sets up initial conditions.
Array< OneD, Array< OneD, NekDouble > > m_beta_trace
Array< OneD, Array< OneD, NekDouble > > m_W2
Array< OneD, Array< OneD, NekDouble > > m_trace_fwd_normal
PulseWavePressureAreaSharedPtr m_pressureArea
void CheckPoint_Output(const int n)
Array< OneD, Array< OneD, NekDouble > > m_W1
virtual void v_InitObject(bool DeclareField=false) override
std::vector< int > m_domOrder
void MergingRiemann(Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0, Array< OneD, NekDouble > &alpha)
Riemann Problem for Merging Flow.
void InterfaceRiemann(Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0, Array< OneD, NekDouble > &alpha)
Riemann Problem for Interface/Junction.
Array< OneD, Array< OneD, NekDouble > > m_alpha_trace
std::map< int, SpatialDomains::CompositeMap > m_domain
UpwindTypePulse m_upwindTypePulse
Array< OneD, Array< OneD, NekDouble > > m_A_0_trace
Array< OneD, Array< OneD, NekDouble > > m_pressure
virtual ~PulseWaveSystem()
Destructor.
virtual NekDouble v_LinfError(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray) override
Compute the L_inf error between fields and a given exact solution.
Array< OneD, Array< OneD, NekDouble > > m_gamma
void WriteVessels(const std::string &outname)
Write input fields to the given filename.
void GetCommArray(std::map< int, LibUtilities::CommSharedPtr > &retval)
Set and retrn a series of communicators for each partition.
Array< OneD, Array< OneD, NekDouble > > m_alpha
Gs::gs_data * m_intComm
Communicator for interfaces.
std::vector< std::vector< InterfacePointShPtr > > m_vesselIntfcs
Array< OneD, Array< OneD, NekDouble > > m_PWV
virtual void v_Output(void) override
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
void FillDataFromInterfacePoint(InterfacePointShPtr &I, const Array< OneD, const Array< OneD, NekDouble > > &field, NekDouble &A, NekDouble &u, NekDouble &beta, NekDouble &A_0, NekDouble &alpha)
std::vector< std::vector< InterfacePointShPtr > > m_bifurcations
void SetUpDomainInterfaceBCs(SpatialDomains::BoundaryConditions &AllBcs, std::map< int, LibUtilities::CommSharedPtr > &domComm)
std::vector< std::vector< InterfacePointShPtr > > m_mergingJcts
virtual NekDouble v_L2Error(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false) override
Compute the L2 error between fields and a given exact solution.
virtual void v_DoSolve() override
Solves an unsteady problem.
void EnforceInterfaceConditions(const Array< OneD, const Array< OneD, NekDouble > > &fields)
Array< OneD, Array< OneD, NekDouble > > m_beta
PulseWaveSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises PulseWaveSystem class members.
void BifurcationRiemann(Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0, Array< OneD, NekDouble > &alpha)
Riemann Problem for Bifurcation.
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_time
Current time of simulation.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void EvaluateExactSolution(int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
Evaluates an exact solution.
std::string m_sessionName
Name of the session.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
int m_steps
Number of steps to take.
int m_NumQuadPointsError
Number of Quadrature points used to work out the error.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
void AddBoundaryRegions(const int regionID, BoundaryRegionShPtr &bRegion)
Definition: Conditions.h:239
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:234
void AddBoundaryConditions(const int regionID, BoundaryConditionMapShPtr &bCond)
Definition: Conditions.h:249
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:270
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:192
@ gs_add
Definition: GsLib.hpp:62
@ gs_max
Definition: GsLib.hpp:65
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::map< int, CompositeSharedPtr > BoundaryRegion
Definition: Conditions.h:207
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
std::shared_ptr< BoundaryRegion > BoundaryRegionShPtr
Definition: Conditions.h:208
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::map< std::string, BoundaryConditionShPtr > BoundaryConditionMap
Definition: Conditions.h:217
std::shared_ptr< Composite > CompositeSharedPtr
Definition: MeshGraph.h:137
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::shared_ptr< BoundaryConditionMap > BoundaryConditionMapShPtr
Definition: Conditions.h:218
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:54
std::vector< double > d(NPUPPER *NPUPPER)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< InterfacePoint > InterfacePointShPtr
PressureAreaFactory & GetPressureAreaFactory()
void Multiply(NekMatrix< ResultDataType, StandardMatrixTag > &result, const NekMatrix< LhsDataType, LhsMatrixType > &lhs, const ResultDataType &rhs)
@ SIZE_UpwindTypePulse
Length of enum list.
const char *const UpwindTypeMapPulse[]
double NekDouble
DataType Dot(const NekVector< DataType > &lhs, const NekVector< DataType > &rhs)
Definition: NekVector.cpp:1193
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:890
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:1018
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.cpp:940
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294