Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DisContField1D.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 //
3 // File DisContField1D.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: Field definition for 1D domain with boundary
33 // conditions using LDG-H
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 #include <StdRegions/StdSegExp.h>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46  namespace MultiRegions
47  {
48  /**
49  * @class DisContField1D
50  * This class augments the list of local expansions inherited from
51  * ExpList1D with boundary conditions. Inter-element boundaries are
52  * handled using an discontinuous Galerkin scheme.
53  */
54 
55  /**
56  * Constructs an empty expansion list with no boundary conditions.
57  */
58  DisContField1D::DisContField1D():
59  ExpList1D(),
60  m_bndCondExpansions(),
61  m_bndConditions()
62  {
63  }
64 
65 
66  /**
67  * An expansion list for the boundary expansions is generated first for
68  * the field. These are subsequently evaluated for time zero. The trace
69  * map is then constructed.
70  * @param graph1D A mesh, containing information about the domain
71  * and the spectral/hp element expansions.
72  * @param bcs Information about the enforced boundary
73  * conditions.
74  * @param variable The session variable associated with the
75  * boundary conditions to enforce.
76  * @param solnType Type of global system to use.
77  */
81  const std::string &variable,
82  const bool SetUpJustDG)
83  : ExpList1D(pSession,graph1D),
84  m_bndCondExpansions(),
85  m_bndConditions()
86  {
88 
89  GenerateBoundaryConditionExpansion(graph1D,bcs,variable);
90  EvaluateBoundaryConditions(0.0, variable);
91  ApplyGeomInfo();
92  FindPeriodicVertices(bcs,variable);
93 
94  if(SetUpJustDG)
95  {
96  SetUpDG(variable);
97  }
98  else
99  {
100  int i;
101  Array<OneD, int> ElmtID, VertexID;
102  GetBoundaryToElmtMap(ElmtID, VertexID);
103 
104  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
105  {
107  locExpList = m_bndCondExpansions[i];
108 
110  (*m_exp)[ElmtID[i]]->
111  as<LocalRegions::Expansion1D>();
113  locExpList->GetExp(0)->
114  as<LocalRegions::Expansion0D>();
115 
116  exp0d->SetAdjacentElementExp(VertexID[i], exp1d);
117  }
118 
120  }
121 
122  }
123 
124  void DisContField1D::SetUpDG(const std::string &variable)
125  {
126  // Check for multiple calls
128  {
129  return;
130  }
131 
133  boost::dynamic_pointer_cast<SpatialDomains::MeshGraph1D>(
134  m_graph);
135 
137 
142  *m_exp,graph1D,
144 
145  m_trace = boost::dynamic_pointer_cast<ExpList>(trace);
146 
148  AllocateSharedPtr(m_session, graph1D, trace, *this,
150  m_periodicVerts, variable);
151 
152  // Scatter trace points to 1D elements. For each element, we find
153  // the trace point associated to each vertex. The element then
154  // retains a pointer to the trace space points, to ensure
155  // uniqueness of normals when retrieving from two adjoining
156  // elements which do not lie in a plane.
157 
158  int ElmtPointGeom = 0;
159  int TracePointGeom = 0;
162  for (int i = 0; i < m_exp->size(); ++i)
163  {
164  exp1d = (*m_exp)[i]->as<LocalRegions::Expansion1D>();
165  for (int j = 0; j < exp1d->GetNverts(); ++j)
166  {
167  ElmtPointGeom = (exp1d->GetGeom1D())->GetVid(j);
168 
169  for (int k = 0; k < m_trace->GetExpSize(); ++k)
170  {
171  TracePointGeom = m_trace->GetExp(k)->GetGeom()->GetVid(0);
172 
173  if (TracePointGeom == ElmtPointGeom)
174  {
175  exp0d = m_trace->GetExp(k)->as<LocalRegions::Expansion0D>();
176  exp0d->SetAdjacentElementExp(j,exp1d);
177  break;
178  }
179  }
180  }
181  }
182 
184 
185  int cnt, n, e;
186 
187  // Identify boundary verts
188  for(cnt = 0, n = 0; n < m_bndCondExpansions.num_elements(); ++n)
189  {
190  if (m_bndConditions[n]->GetBoundaryConditionType() !=
192  {
193  for(e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
194  {
195  m_boundaryVerts.insert(
196  m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
197  }
198  }
199  else
200  {
201  ASSERTL0(false,"Periodic verts need setting up");
202  }
203  cnt += m_bndCondExpansions[n]->GetExpSize();
204  }
205 
206  // Set up left-adjacent edge list.
207  m_leftAdjacentVerts.resize(2*((*m_exp).size()));
208 
209  // count size of trace
210  for (cnt = n = 0; n < m_exp->size(); ++n)
211  {
212  for (int v = 0; v < (*m_exp)[n]->GetNverts(); ++v, ++cnt)
213  {
215  }
216  }
217 
218 
219  boost::unordered_map<int,pair<int,int> > perVertToExpMap;
220  boost::unordered_map<int,pair<int,int> >::iterator it2;
221  for (n = 0; n < m_exp->size(); ++n)
222  {
223  for (int v = 0; v < (*m_exp)[n]->GetNverts(); ++v)
224  {
226  (*m_exp)[n]->GetGeom()->GetVid(v));
227 
228  if (it != m_periodicVerts.end())
229  {
230  perVertToExpMap[it->first] = make_pair(n,v);
231  }
232  }
233  }
234 
235 
236  // Set up mapping to copy Fwd of periodic bcs to Bwd of other edge.
237  for (n = 0; n < m_exp->size(); ++n)
238  {
239  for (int v = 0; v < (*m_exp)[n]->GetNverts(); ++v)
240  {
241  int vertGeomId = (*m_exp)[n]->GetGeom()->GetVid(v);
242 
243  // Check to see if this face is periodic.
244  PeriodicMap::iterator it = m_periodicVerts.find(vertGeomId);
245 
246  if (it != m_periodicVerts.end())
247  {
248  const PeriodicEntity &ent = it->second[0];
249  it2 = perVertToExpMap.find(ent.id);
250 
251  if (it2 == perVertToExpMap.end())
252  {
253  if (m_session->GetComm()->GetRowComm()->GetSize() > 1 &&
254  !ent.isLocal)
255  {
256  continue;
257  }
258  else
259  {
260  ASSERTL1(false, "Periodic vert not found!");
261  }
262  }
263 
264  int offset = m_trace->GetPhys_Offset((m_traceMap->GetElmtToTrace())
265  [n][v]->GetElmtId());
266  m_periodicFwdCopy.push_back(offset);
267 
268  int offset2 = m_trace->GetPhys_Offset((m_traceMap->GetElmtToTrace())
269  [it2->second.first]
270  [it2->second.second]->GetElmtId());
271  m_periodicBwdCopy.push_back(offset2);
272  }
273  }
274  }
275  }
276 
277 
278  bool DisContField1D::IsLeftAdjacentVertex(const int n, const int e)
279  {
282  m_traceMap->GetElmtToTrace()[n][e]->as<LocalRegions::Expansion0D>();
283 
284 
285  bool fwd = true;
286  if (traceEl->GetLeftAdjacentElementVertex () == -1 ||
287  traceEl->GetRightAdjacentElementVertex() == -1)
288  {
289  // Boundary edge (1 connected element). Do nothing in
290  // serial.
291  it = m_boundaryVerts.find(traceEl->GetElmtId());
292 
293  // If the edge does not have a boundary condition set on
294  // it, then assume it is a partition edge or periodic.
295  if (it == m_boundaryVerts.end())
296  {
297  int traceGeomId = traceEl->GetGeom0D()->GetGlobalID();
299  traceGeomId);
300 
301  if (pIt != m_periodicVerts.end() && !pIt->second[0].isLocal)
302  {
303  fwd = traceGeomId == min(traceGeomId,pIt->second[0].id);
304  }
305  else
306  {
307  int offset = m_trace->GetPhys_Offset(traceEl->GetElmtId());
308  fwd = m_traceMap->
309  GetTraceToUniversalMapUnique(offset) >= 0;
310  }
311  }
312  }
313  else if (traceEl->GetLeftAdjacentElementVertex () != -1 &&
314  traceEl->GetRightAdjacentElementVertex() != -1)
315  {
316  // Non-boundary edge (2 connected elements).
317  fwd = dynamic_cast<Nektar::StdRegions::StdExpansion*>
318  (traceEl->GetLeftAdjacentElementExp().get()) ==
319  (*m_exp)[n].get();
320  }
321  else
322  {
323  ASSERTL2(false, "Unconnected trace element!");
324  }
325 
326  return fwd;
327  }
328 
329 
330  // Given all boundary regions for the whole solution determine
331  // which ones (if any) are part of domain and ensure all other
332  // conditions are given as UserDefined Dirichlet.
334  const SpatialDomains::CompositeMap &domain,
336  const std::string &variable)
337  {
339 
341 
343  map<int,int> GeometryToRegionsMap;
344 
345  SpatialDomains::BoundaryRegionCollection::const_iterator it;
346 
348  = Allbcs.GetBoundaryRegions();
350  = Allbcs.GetBoundaryConditions();
351 
352  // Set up a map of all boundary regions
353  for(it = bregions.begin(); it != bregions.end(); ++it)
354  {
356  for (bregionIt = it->second->begin();
357  bregionIt != it->second->end(); bregionIt++)
358  {
359  // can assume that all regions only contain one point in 1D
360  // Really do not need loop above
361  int id = (*(bregionIt->second))[0]->GetGlobalID();
362  GeometryToRegionsMap[id] = it->first;
363  }
364  }
365 
367  map<int,SpatialDomains::GeometrySharedPtr> EndOfDomain;
368 
369  // Now find out which points in domain have only one vertex
370  for(domIt = domain.begin(); domIt != domain.end(); ++domIt)
371  {
372  SpatialDomains::Composite geomvector = domIt->second;
373  for(int i = 0; i < geomvector->size(); ++i)
374  {
375  for(int j = 0; j < 2; ++j)
376  {
377  int vid = (*geomvector)[i]->GetVid(j);
378  if(EndOfDomain.count(vid) == 0)
379  {
380  EndOfDomain[vid] = (*geomvector)[i]->GetVertex(j);
381  }
382  else
383  {
384  EndOfDomain.erase(vid);
385  }
386  }
387  }
388  }
389  ASSERTL1(EndOfDomain.size() == 2,"Did not find two ends of domain");
390 
392  int numNewBc = 1;
393  for(regIt = EndOfDomain.begin(); regIt != EndOfDomain.end(); ++regIt)
394  {
395  if(GeometryToRegionsMap.count(regIt->first) != 0) // Set up boundary condition up
396  {
397  map<int,int>::iterator iter = GeometryToRegionsMap.find(regIt->first);
398  ASSERTL1(iter != GeometryToRegionsMap.end(),"Failied to find GeometryToRegionMap");
399  int regionId = iter->second;
400  SpatialDomains::BoundaryRegionCollection::const_iterator bregionsIter = bregions.find(regionId);
401  ASSERTL1(bregionsIter != bregions.end(),"Failed to find boundary region");
402  SpatialDomains::BoundaryRegionShPtr breg = bregionsIter->second;
403  returnval->AddBoundaryRegions (regionId,breg);
404 
405  SpatialDomains::BoundaryConditionCollection::const_iterator bconditionsIter = bconditions.find(regionId);
406  ASSERTL1(bconditionsIter != bconditions.end(),"Failed to find boundary collection");
407  SpatialDomains::BoundaryConditionMapShPtr bcond = bconditionsIter->second;
408  returnval->AddBoundaryConditions(regionId,bcond);
409  }
410  else // Set up an undefined region.
411  {
413 
414  // Set up Composite (GemetryVector) to contain vertex and put into bRegion
416  gvec->push_back(regIt->second);
417  (*breg)[regIt->first] = gvec;
418 
419  returnval->AddBoundaryRegions(bregions.size()+numNewBc,breg);
420 
422 
423  // Set up just boundary condition for this variable.
425  (*bCondition)[variable] = notDefinedCondition;
426 
427  returnval->AddBoundaryConditions(bregions.size()+numNewBc,bCondition);
428  ++numNewBc;
429 
430  }
431  }
432 
433  return returnval;
434  }
435 
436  /**
437  * Constructor for use in multidomain computations where a
438  * domain list can be passed instead of graph1D
439  *
440  * @param domain Subdomain specified in the inputfile from
441  * which the DisContField1D is set up
442  */
444  const LibUtilities::SessionReaderSharedPtr &pSession,
445  const SpatialDomains::MeshGraphSharedPtr &graph1D,
446  const SpatialDomains::CompositeMap &domain,
447  const SpatialDomains::BoundaryConditions &Allbcs,
448  const std::string &variable,
449  bool SetToOneSpaceDimension):
450  ExpList1D(pSession,graph1D,domain, true,variable,SetToOneSpaceDimension),
451  m_bndCondExpansions(),
452  m_bndConditions()
453  {
454  SpatialDomains::BoundaryConditionsSharedPtr DomBCs = GetDomainBCs(domain,Allbcs,variable);
455 
456  GenerateBoundaryConditionExpansion(graph1D,*DomBCs,variable);
457  EvaluateBoundaryConditions(0.0, variable);
458  ApplyGeomInfo();
459  FindPeriodicVertices(*DomBCs,variable);
460 
461  SetUpDG(variable);
462  }
463 
464  /**
465  * Constructs a field as a copy of an existing field.
466  * @param In Existing DisContField1D object to copy.
467  */
469  ExpList1D(In),
470  m_bndCondExpansions(In.m_bndCondExpansions),
471  m_bndConditions(In.m_bndConditions),
472  m_globalBndMat(In.m_globalBndMat),
473  m_trace(In.m_trace),
474  m_traceMap(In.m_traceMap),
475  m_boundaryVerts(In.m_boundaryVerts),
476  m_periodicVerts(In.m_periodicVerts),
477  m_periodicFwdCopy(In.m_periodicFwdCopy),
478  m_periodicBwdCopy(In.m_periodicBwdCopy),
479  m_leftAdjacentVerts(In.m_leftAdjacentVerts)
480  {
481  }
482 
483 
484  /**
485  * Constructs a field as a copy of an existing explist1D field.
486  * @param In Existing ExpList1D object to copy.
487  */
489  ExpList1D(In)
490  {
491  }
492 
493  /**
494  *
495  */
497  {
498  }
499 
500 
501  /**
502  * Generate the boundary condition expansion list
503  * @param graph1D A mesh, containing information about the domain
504  * and the spectral/hp element expansions.
505  * @param bcs Information about the enforced boundary
506  * conditions.
507  * @param variable The session variable associated with the
508  * boundary conditions to enforce.
509  */
511  const SpatialDomains::MeshGraphSharedPtr &graph1D,
513  const std::string variable)
514  {
515  int cnt = 0;
516 
518  = bcs.GetBoundaryRegions();
520  = bcs.GetBoundaryConditions();
521  SpatialDomains::BoundaryRegionCollection::const_iterator it;
522 
523  // count the number of non-periodic boundary points
524  for (it = bregions.begin(); it != bregions.end(); ++it)
525  {
526  const SpatialDomains::BoundaryConditionShPtr boundaryCondition =
527  GetBoundaryCondition(bconditions, it->first, variable);
528  if (boundaryCondition->GetBoundaryConditionType() !=
530  {
532  for (bregionIt = it->second->begin();
533  bregionIt != it->second->end(); bregionIt++)
534  {
535  cnt += bregionIt->second->size();
536  }
537  }
538  }
539 
542 
545 
546  SetBoundaryConditionExpansion(graph1D,bcs,variable,
548  m_bndConditions);
549  }
550 
551 
552  /**
553  * @param graph1D A mesh containing information about the domain
554  * and the spectral/hp element expansion.
555  * @param bcs Information about the boundary conditions.
556  * @param variable Specifies the field.
557  * @param periodicVertices Map into which the list of periodic
558  * vertices is placed.
559  */
562  const std::string variable)
563  {
565  = bcs.GetBoundaryRegions();
567  = bcs.GetBoundaryConditions();
568 
570  = boost::dynamic_pointer_cast<
572  SpatialDomains::BoundaryRegionCollection::const_iterator it;
573 
575  m_session->GetComm()->GetRowComm();
576 
577  int i, region1ID, region2ID;
578 
580 
581  map<int,int> BregionToVertMap;
582 
583  // Construct list of all periodic Region and their global vertex on
584  // this process.
585  for (it = bregions.begin(); it != bregions.end(); ++it)
586  {
587  locBCond = GetBoundaryCondition(bconditions, it->first, variable);
588 
589  if (locBCond->GetBoundaryConditionType()
591  {
592  continue;
593  }
594 
595  int id = (*(it->second->begin()->second))[0]->GetGlobalID();
596 
597  BregionToVertMap[it->first] = id;
598  }
599 
601  set<int> islocal;
602 
603  int n = vComm->GetSize();
604  int p = vComm->GetRank();
605 
606  Array<OneD, int> nregions(n, 0);
607  nregions[p] = BregionToVertMap.size();
608  vComm->AllReduce(nregions, LibUtilities::ReduceSum);
609 
610  int totRegions = Vmath::Vsum(n, nregions, 1);
611 
612  Array<OneD, int> regOffset(n, 0);
613 
614  for (i = 1; i < n; ++i)
615  {
616  regOffset[i] = regOffset[i-1] + nregions[i-1];
617  }
618 
619  Array<OneD, int> bregmap(totRegions, 0);
620  Array<OneD, int> bregid (totRegions, 0);
621  for(i = regOffset[p], iit = BregionToVertMap.begin();
622  iit != BregionToVertMap.end(); ++iit, ++i)
623  {
624  bregid [i] = iit->first;
625  bregmap[i] = iit->second;
626  islocal.insert(iit->first);
627  }
628 
629  vComm->AllReduce(bregmap, LibUtilities::ReduceSum);
630  vComm->AllReduce(bregid, LibUtilities::ReduceSum);
631 
632  for (int i = 0; i < totRegions; ++i)
633  {
634  BregionToVertMap[bregid[i]] = bregmap[i];
635  }
636 
637  // Construct list of all periodic pairs local to this process.
638  for (it = bregions.begin(); it != bregions.end(); ++it)
639  {
640  locBCond = GetBoundaryCondition(bconditions, it->first, variable);
641 
642  if (locBCond->GetBoundaryConditionType()
644  {
645  continue;
646  }
647 
648  // Identify periodic boundary region IDs.
649  region1ID = it->first;
650  region2ID = boost::static_pointer_cast<
652  locBCond)->m_connectedBoundaryRegion;
653 
654  ASSERTL0(BregionToVertMap.count(region1ID) != 0,
655  "Cannot determine vertex of region1ID");
656 
657  ASSERTL0(BregionToVertMap.count(region2ID) != 0,
658  "Cannot determine vertex of region2ID");
659 
660  PeriodicEntity ent(BregionToVertMap[region2ID],
662  islocal.count(region2ID) != 0);
663 
664  m_periodicVerts[BregionToVertMap[region1ID]].push_back(ent);
665  }
666  }
667 
668 
669  /**
670  * @param graph1D A mesh containing information about the domain
671  * and the Spectral/hp element expansion.
672  * @param bcs Information about the boundary conditions.
673  * @param variable Specifies the field.
674  * @param bndCondExpansions Array of ExpList1D objects each
675  * containing a 1D spectral/hp element expansion
676  * on a single boundary region.
677  * @param bncConditions Array of BoundaryCondition objects which
678  * contain information about the boundary
679  * conditions on the different boundary regions.
680  */
682  const SpatialDomains::MeshGraphSharedPtr &graph1D,
684  const std::string variable,
686  &bndCondExpansions,
687  Array<OneD, SpatialDomains
688  ::BoundaryConditionShPtr> &bndConditions)
689  {
690  int k;
691  int cnt = 0;
692 
694  = bcs.GetBoundaryRegions();
696  = bcs.GetBoundaryConditions();
697  SpatialDomains::BoundaryRegionCollection::const_iterator it;
698 
702 
703  cnt = 0;
704  // list Dirichlet boundaries first
705  for (it = bregions.begin(); it != bregions.end(); ++it)
706  {
707  locBCond = GetBoundaryCondition(
708  bconditions, it->first, variable);
709 
710  if (locBCond->GetBoundaryConditionType() ==
712 
713  {
715  for (bregionIt = it->second->begin();
716  bregionIt != it->second->end(); bregionIt++)
717  {
718  for (k = 0; k < bregionIt->second->size(); k++)
719  {
720  if ((vert = boost::dynamic_pointer_cast
722  (*bregionIt->second)[k])))
723  {
724  locPointExp
727  bndCondExpansions[cnt] = locPointExp;
728  bndConditions[cnt++] = locBCond;
729  }
730  else
731  {
732  ASSERTL0(false,
733  "dynamic cast to a vertex failed");
734  }
735  }
736  }
737  }
738  } // end if Dirichlet
739 
740  // then, list the other (non-periodic) boundaries
741  for (it = bregions.begin(); it != bregions.end(); ++it)
742  {
743  locBCond = GetBoundaryCondition(bconditions, it->first, variable);
744 
745  switch(locBCond->GetBoundaryConditionType())
746  {
749  case SpatialDomains::eNotDefined: // presume this will be reused as Neuman, Robin or Dirichlet later
750  {
752  for (bregionIt = it->second->begin();
753  bregionIt != it->second->end(); bregionIt++)
754  {
755  for (k = 0; k < bregionIt->second->size(); k++)
756  {
757  if((vert = boost::dynamic_pointer_cast
759  (*bregionIt->second)[k])))
760  {
761  locPointExp
764  bndCondExpansions[cnt] = locPointExp;
765  bndConditions[cnt++] = locBCond;
766  }
767  else
768  {
769  ASSERTL0(false,
770  "dynamic cast to a vertex failed");
771  }
772  }
773  }
774  }
775  // do nothing for these types
778  break;
779  default:
780  ASSERTL0(false,"This type of BC not implemented yet");
781  break;
782  }
783  }
784  }
785 
786  /**
787  *
788  */
790  const GlobalLinSysKey &mkey)
791  {
793  "Routine currently only tested for HybridDGHelmholtz");
794 
796  "Full matrix global systems are not supported for HDG "
797  "expansions");
798 
800  ==m_traceMap->GetGlobalSysSolnType(),
801  "The local to global map is not set up for the requested "
802  "solution type");
803 
804  GlobalLinSysSharedPtr glo_matrix;
805  GlobalLinSysMap::iterator matrixIter = m_globalBndMat->find(mkey);
806 
807  if (matrixIter == m_globalBndMat->end())
808  {
809  glo_matrix = GenGlobalBndLinSys(mkey,m_traceMap);
810  (*m_globalBndMat)[mkey] = glo_matrix;
811  }
812  else
813  {
814  glo_matrix = matrixIter->second;
815  }
816 
817  return glo_matrix;
818  }
819 
820 
822  {
823  if(m_negatedFluxNormal.size() == 0)
824  {
826  &elmtToTrace = m_traceMap->GetElmtToTrace();
827 
828  m_negatedFluxNormal.resize(2*GetExpSize());
829 
830  for(int i = 0; i < GetExpSize(); ++i)
831  {
832 
833  for(int v = 0; v < 2; ++v)
834  {
835 
837  elmtToTrace[i][v]->as<LocalRegions::Expansion0D>();
838 
839  if(vertExp->GetLeftAdjacentElementExp()->GetGeom()->GetGlobalID() != (*m_exp)[i]->GetGeom()->GetGlobalID())
840  {
841  m_negatedFluxNormal[2*i+v] = true;
842  }
843  else
844  {
845  m_negatedFluxNormal[2*i+v] = false;
846  }
847  }
848  }
849 
850  }
851 
852  return m_negatedFluxNormal;
853  }
854 
855  /**
856  * Generate the forward or backward state for each trace point.
857  * @param Fwd Forward state.
858  * @param Bwd Backward state.
859  */
862  {
863  v_GetFwdBwdTracePhys(m_phys,Fwd,Bwd);
864  }
865 
866 
867  /**
868  * @brief This method extracts the "forward" and "backward" trace data
869  * from the array @a field and puts the data into output vectors @a Fwd
870  * and @a Bwd.
871  *
872  * We first define the convention which defines "forwards" and
873  * "backwards". First an association is made between the vertex of each
874  * element and its corresponding vertex in the trace space using the
875  * mapping #m_traceMap. The element can either be left-adjacent or
876  * right-adjacent to this trace edge (see
877  * Expansion0D::GetLeftAdjacentElementExp). Boundary edges are never
878  * left-adjacent since elemental left-adjacency is populated first.
879  *
880  * If the element is left-adjacent we extract the vertex trace data from
881  * @a field into the forward trace space @a Fwd; otherwise, we place it
882  * in the backwards trace space @a Bwd. In this way, we form a unique
883  * set of trace normals since these are always extracted from
884  * left-adjacent elements.
885  *
886  * @param field is a NekDouble array which contains the 1D data
887  * from which we wish to extract the backward and forward
888  * orientated trace/edge arrays.
889  * @param Fwd The resulting forwards space.
890  * @param Bwd The resulting backwards space.
891  */
892 
894  const Array<OneD, const NekDouble> &field,
897  {
898  // Counter variables
899  int n, v;
900 
901  // Number of elements
902  int nElements = GetExpSize();
903 
904  // Initial index of each element
905  int phys_offset;
906 
908  &elmtToTrace = m_traceMap->GetElmtToTrace();
909 
910  // Set forward and backard state to zero
911  Vmath::Zero(Fwd.num_elements(), Fwd, 1);
912  Vmath::Zero(Bwd.num_elements(), Bwd, 1);
913 
914  int cnt;
915 
916  // Loop on the elements
917  for (cnt = n = 0; n < nElements; ++n)
918  {
919  // Set the offset of each element
920  phys_offset = GetPhys_Offset(n);
921 
922  for(v = 0; v < 2; ++v, ++cnt)
923  {
924  int offset = m_trace->GetPhys_Offset(elmtToTrace[n][v]->GetElmtId());
925 
926  if (m_leftAdjacentVerts[cnt])
927  {
928  (*m_exp)[n]->GetVertexPhysVals(v, field + phys_offset,
929  Fwd[offset]);
930  }
931  else
932  {
933  (*m_exp)[n]->GetVertexPhysVals(v, field + phys_offset,
934  Bwd[offset]);
935  }
936  }
937  }
938 
939  // Fill boundary conditions into missing elements.
940  int id = 0;
941 
942  for(cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
943  {
944  if (m_bndConditions[n]->GetBoundaryConditionType() ==
946  {
947  id = m_trace->GetPhys_Offset(m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt));
948  Bwd[id] = m_bndCondExpansions[n]->GetPhys()[0]; //this is not getting the correct value?
949  cnt++;
950  }
951  else if (m_bndConditions[n]->GetBoundaryConditionType() ==
953  m_bndConditions[n]->GetBoundaryConditionType() ==
955  {
956  ASSERTL0((m_bndCondExpansions[n]->GetPhys())[0]==0.0,
957  "Method not set up for non-zero Neumann "
958  "boundary condition");
959  id = m_trace->GetPhys_Offset(m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt));
960  Bwd[id] = Fwd[id];
961 
962  cnt++;
963  }
964  else if (m_bndConditions[n]->GetBoundaryConditionType() ==
966  {
967  // Do nothing
968  }
969  else if (m_bndConditions[n]->GetBoundaryConditionType() !=
971  {
972  ASSERTL0(false,
973  "Method not set up for this boundary condition.");
974  }
975  }
976 
977  // Copy any periodic boundary conditions.
978  for (n = 0; n < m_periodicFwdCopy.size(); ++n)
979  {
980  Bwd[m_periodicBwdCopy[n]] = Fwd[m_periodicFwdCopy[n]];
981  }
982 
983  // Do parallel exchange for forwards/backwards spaces.
984  m_traceMap->UniversalTraceAssemble(Fwd);
985  m_traceMap->UniversalTraceAssemble(Bwd);
986 
987  }
988 
989 
991  Array<OneD, NekDouble> &outarray)
992  {
993  ASSERTL1(m_physState == true,"local physical space is not true ");
994  v_ExtractTracePhys(m_phys, outarray);
995  }
996 
997  /**
998  * @brief This method extracts the trace (verts in 1D) from the field @a
999  * inarray and puts the values in @a outarray.
1000  *
1001  * It assumes the field is C0 continuous so that it can overwrite the
1002  * edge data when visited by the two adjacent elements.
1003  *
1004  * @param inarray An array containing the 1D data from which we wish
1005  * to extract the edge data.
1006  * @param outarray The resulting edge information.
1007  *
1008  * This will not work for non-boundary expansions
1009  */
1011  const Array<OneD, const NekDouble> &inarray,
1012  Array<OneD, NekDouble> &outarray)
1013  {
1014  // Loop over elemente and collect forward expansion
1015  int nexp = GetExpSize();
1016  int n,p,offset,phys_offset;
1017 
1018  ASSERTL1(outarray.num_elements() >= m_trace->GetExpSize(),
1019  "input array is of insufficient length");
1020 
1021  for (n = 0; n < nexp; ++n)
1022  {
1023  phys_offset = GetPhys_Offset(n);
1024 
1025  for (p = 0; p < (*m_exp)[n]->GetNverts(); ++p)
1026  {
1027  offset = m_trace->GetPhys_Offset(
1028  (m_traceMap->GetElmtToTrace())[n][p]->GetElmtId());
1029  (*m_exp)[n]->GetVertexPhysVals(p, inarray + phys_offset,
1030  outarray[offset]);
1031  }
1032  }
1033  }
1034 
1036  const Array<OneD, const NekDouble> &Fn,
1037  Array<OneD, NekDouble> &outarray)
1038  {
1039  int n,offset, t_offset;
1040 
1042  &elmtToTrace = m_traceMap->GetElmtToTrace();
1043 
1044  vector<bool> negatedFluxNormal = GetNegatedFluxNormal();
1045 
1046  for (n = 0; n < GetExpSize(); ++n)
1047  {
1048  // Number of coefficients on each element
1049  int e_ncoeffs = (*m_exp)[n]->GetNcoeffs();
1050 
1051  offset = GetCoeff_Offset(n);
1052 
1053  // Implementation for every points except Gauss points
1054  if ((*m_exp)[n]->GetBasis(0)->GetBasisType() !=
1056  {
1057  t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][0]->GetElmtId());
1058  if(negatedFluxNormal[2*n])
1059  {
1060  outarray[offset] -= Fn[t_offset];
1061  }
1062  else
1063  {
1064  outarray[offset] += Fn[t_offset];
1065  }
1066 
1067  t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][1]->GetElmtId());
1068 
1069  if(negatedFluxNormal[2*n+1])
1070  {
1071  outarray[offset+(*m_exp)[n]->GetVertexMap(1)] -= Fn[t_offset];
1072  }
1073  else
1074  {
1075  outarray[offset+(*m_exp)[n]->GetVertexMap(1)] += Fn[t_offset];
1076  }
1077 
1078  }
1079  else
1080  {
1081 #if 0
1082  DNekMatSharedPtr m_Ixm;
1085  BS_p(e_ncoeffs,LibUtilities::eGaussGaussLegendre);
1087  BS_k(LibUtilities::eGauss_Lagrange,e_ncoeffs,BS_p);
1088 
1089  BASE = LibUtilities::BasisManager()[BS_k];
1090 
1091  Array<OneD, NekDouble> coords(3, 0.0);
1092 
1093  int j;
1094 
1095  for(p = 0; p < 2; ++p)
1096  {
1097  NekDouble vertnorm = 0.0;
1098  for (int i=0; i<((*m_exp)[n]->
1099  GetVertexNormal(p)).num_elements(); i++)
1100  {
1101  vertnorm += ((*m_exp)[n]->GetVertexNormal(p))[i][0];
1102  coords[0] = vertnorm ;
1103  }
1104 
1105  t_offset = GetTrace()->GetPhys_Offset(n+p);
1106 
1107  if (vertnorm >= 0.0)
1108  {
1109  m_Ixm = BASE->GetI(coords);
1110 
1111 
1112  for (j = 0; j < e_ncoeffs; j++)
1113  {
1114  outarray[offset + j] +=
1115  (m_Ixm->GetPtr())[j] * Fn[t_offset];
1116  }
1117  }
1118 
1119  if (vertnorm < 0.0)
1120  {
1121  m_Ixm = BASE->GetI(coords);
1122 
1123  for (j = 0; j < e_ncoeffs; j++)
1124  {
1125  outarray[offset + j] -=
1126  (m_Ixm->GetPtr())[j] * Fn[t_offset];
1127  }
1128  }
1129  }
1130 #else
1131  int j;
1132  static DNekMatSharedPtr m_Ixm, m_Ixp;
1133  static int sav_ncoeffs = 0;
1134  if(!m_Ixm.get() || e_ncoeffs != sav_ncoeffs)
1135  {
1138  BS_p(e_ncoeffs,LibUtilities::eGaussGaussLegendre);
1140  BS_k(LibUtilities::eGauss_Lagrange,e_ncoeffs,BS_p);
1141 
1142  BASE = LibUtilities::BasisManager()[BS_k];
1143 
1144  Array<OneD, NekDouble> coords(1, 0.0);
1145 
1146  coords[0] = -1.0;
1147  m_Ixm = BASE->GetI(coords);
1148 
1149  coords[0] = 1.0;
1150  m_Ixp = BASE->GetI(coords);
1151 
1152  sav_ncoeffs = e_ncoeffs;
1153  }
1154 
1155  t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][0]->GetElmtId());
1156  if(negatedFluxNormal[2*n])
1157  {
1158  for (j = 0; j < e_ncoeffs; j++)
1159  {
1160  outarray[offset + j] -=
1161  (m_Ixm->GetPtr())[j] * Fn[t_offset];
1162  }
1163  }
1164  else
1165  {
1166  for (j = 0; j < e_ncoeffs; j++)
1167  {
1168  outarray[offset + j] +=
1169  (m_Ixm->GetPtr())[j] * Fn[t_offset];
1170  }
1171  }
1172 
1173  t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][1]->GetElmtId());
1174  if (negatedFluxNormal[2*n+1])
1175  {
1176  for (j = 0; j < e_ncoeffs; j++)
1177  {
1178  outarray[offset + j] -=
1179  (m_Ixp->GetPtr())[j] * Fn[t_offset];
1180  }
1181  }
1182  else
1183  {
1184  for (j = 0; j < e_ncoeffs; j++)
1185  {
1186  outarray[offset + j] +=
1187  (m_Ixp->GetPtr())[j] * Fn[t_offset];
1188  }
1189  }
1190 #endif
1191  }
1192  }
1193  }
1194 
1195 
1197  const Array<OneD, const NekDouble> &inarray,
1198  Array<OneD, NekDouble> &outarray,
1199  const FlagList &flags,
1200  const StdRegions::ConstFactorMap &factors,
1201  const StdRegions::VarCoeffMap &varcoeff,
1202  const Array<OneD, const NekDouble> &dirForcing)
1203  {
1204  int i,n,cnt,nbndry;
1205  int nexp = GetExpSize();
1207  DNekVec F(m_ncoeffs,f,eWrapper);
1208  Array<OneD,NekDouble> e_f, e_l;
1209 
1210  //----------------------------------
1211  // Setup RHS Inner product
1212  //----------------------------------
1213  IProductWRTBase(inarray,f);
1214  Vmath::Neg(m_ncoeffs,f,1);
1215 
1216  //----------------------------------
1217  // Solve continuous Boundary System
1218  //----------------------------------
1219  int GloBndDofs = m_traceMap->GetNumGlobalBndCoeffs();
1220  int NumDirBCs = m_traceMap->GetNumLocalDirBndCoeffs();
1221  int e_ncoeffs,id;
1222 
1223  GlobalMatrixKey HDGLamToUKey(
1226  factors,
1227  varcoeff);
1228 
1229  const DNekScalBlkMatSharedPtr &HDGLamToU =
1230  GetBlockMatrix(HDGLamToUKey);
1231 
1232  // Retrieve global trace space storage, \Lambda, from trace expansion
1234  (m_traceMap->GetNumLocalBndCoeffs());
1235 
1236 
1237  Array<OneD,NekDouble> BndRhs(GloBndDofs,0.0);
1238  // Zero trace space
1239  Vmath::Zero(GloBndDofs,BndSol,1);
1240 
1241  int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
1242  Array<OneD, NekDouble> loc_lambda(LocBndCoeffs);
1243  DNekVec LocLambda(LocBndCoeffs,loc_lambda,eWrapper);
1244 
1245  //----------------------------------
1246  // Evaluate Trace Forcing
1247  //----------------------------------
1248  // Determing <u_lam,f> terms using HDGLamToU matrix
1249  for (cnt = n = 0; n < nexp; ++n)
1250  {
1251  nbndry = (*m_exp)[n]->NumDGBndryCoeffs();
1252 
1253  e_ncoeffs = (*m_exp)[n]->GetNcoeffs();
1254  e_f = f+m_coeff_offset[n];
1255  e_l = loc_lambda + cnt;
1256 
1257  // use outarray as tmp space
1258  DNekVec Floc (nbndry, e_l, eWrapper);
1259  DNekVec ElmtFce (e_ncoeffs, e_f, eWrapper);
1260  Floc = Transpose(*(HDGLamToU->GetBlock(n,n)))*ElmtFce;
1261 
1262  cnt += nbndry;
1263  }
1264 
1265  // Assemble into global operator
1266  m_traceMap->AssembleBnd(loc_lambda,BndRhs);
1267 
1268  cnt = 0;
1269  // Copy Dirichlet boundary conditions into trace space
1270  for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
1271  {
1272  if (m_bndConditions[i]->GetBoundaryConditionType() ==
1274  {
1275  id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(i);
1276  BndSol[id] = m_bndCondExpansions[i]->GetCoeff(0);
1277  }
1278  else
1279  {
1280  id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(i);
1281  BndRhs[id] += m_bndCondExpansions[i]->GetCoeff(0);
1282  }
1283  }
1284 
1285  //----------------------------------
1286  // Solve trace problem
1287  //----------------------------------
1288  if (GloBndDofs - NumDirBCs > 0)
1289  {
1291  m_traceMap,factors);
1293  LinSys->Solve(BndRhs,BndSol,m_traceMap);
1294  }
1295 
1296  //----------------------------------
1297  // Internal element solves
1298  //----------------------------------
1301  factors);
1302 
1303  const DNekScalBlkMatSharedPtr& InvHDGHelm =
1304  GetBlockMatrix(invHDGhelmkey);
1305  DNekVec out(m_ncoeffs,outarray,eWrapper);
1306  Vmath::Zero(m_ncoeffs,outarray,1);
1307 
1308  // get local trace solution from BndSol
1309  m_traceMap->GlobalToLocalBnd(BndSol,loc_lambda);
1310 
1311  // out = u_f + u_lam = (*InvHDGHelm)*f + (LamtoU)*Lam
1312  out = (*InvHDGHelm)*F + (*HDGLamToU)*LocLambda;
1313  }
1314 
1315  /**
1316  * Based on the expression \f$g(x,t)\f$ for the boundary conditions,
1317  * this function evaluates the boundary conditions for all boundaries
1318  * at time-level \a t.
1319  * @param time The time at which the boundary conditions
1320  * should be evaluated.
1321  * @param bndCondExpansions List of boundary expansions.
1322  * @param bndConditions Information about the boundary conditions.
1323  */
1325  const NekDouble time,
1326  const std::string varName,
1327  const NekDouble x2_in,
1328  const NekDouble x3_in)
1329  {
1330  int i;
1331 
1332  Array<OneD, NekDouble> x0(1);
1333  Array<OneD, NekDouble> x1(1);
1334  Array<OneD, NekDouble> x2(1);
1335 
1336  for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
1337  {
1338  if (time == 0.0 || m_bndConditions[i]->IsTimeDependent())
1339  {
1340  m_bndCondExpansions[i]->GetCoords(x0, x1, x2);
1341 
1342  if (x2_in != NekConstants::kNekUnsetDouble && x3_in !=
1344  {
1345  x1[0] = x2_in;
1346  x2[0] = x3_in;
1347  }
1348 
1349  if (m_bndConditions[i]->GetBoundaryConditionType() ==
1351  {
1352  m_bndCondExpansions[i]->SetCoeff(0,
1353  (boost::static_pointer_cast<SpatialDomains
1354  ::DirichletBoundaryCondition>(m_bndConditions[i])
1355  ->m_dirichletCondition).Evaluate(x0[0],x1[0],x2[0],time));
1356  m_bndCondExpansions[i]->SetPhys(0,m_bndCondExpansions[i]->GetCoeff(0));
1357  }
1358  else if (m_bndConditions[i]->GetBoundaryConditionType()
1360  {
1361  m_bndCondExpansions[i]->SetCoeff(0,
1362  (boost::static_pointer_cast<SpatialDomains
1363  ::NeumannBoundaryCondition>(m_bndConditions[i])
1364  ->m_neumannCondition).Evaluate(x0[0],x1[0],x2[0],time));
1365  }
1366  else if (m_bndConditions[i]->GetBoundaryConditionType()
1368  {
1369  m_bndCondExpansions[i]->SetCoeff(0,
1370  (boost::static_pointer_cast<SpatialDomains
1371  ::RobinBoundaryCondition>(m_bndConditions[i])
1372  ->m_robinFunction).Evaluate(x0[0],x1[0],x2[0],time));
1373 
1374  m_bndCondExpansions[i]->SetPhys(0,
1375  (boost::static_pointer_cast<SpatialDomains
1376  ::RobinBoundaryCondition>(m_bndConditions[i])
1377  ->m_robinPrimitiveCoeff).Evaluate(x0[0],x1[0],x2[0],time));
1378 
1379  }
1380  else if (m_bndConditions[i]->GetBoundaryConditionType()
1382  {
1383  }
1384  else
1385  {
1386  ASSERTL0(false, "This type of BC not implemented yet");
1387  }
1388  }
1389  }
1390  }
1391 
1392  // Set up a list of element ids and edge ids that link to the
1393  // boundary conditions
1395  Array<OneD, int> &ElmtID, Array<OneD,int> &VertID)
1396  {
1397  map<int, int> VertGID;
1398  int i,n,id;
1399  int bid,cnt,Vid;
1400  int nbcs = m_bndConditions.num_elements();
1401 
1402  // make sure arrays are of sufficient length
1403  if (ElmtID.num_elements() != nbcs)
1404  {
1405  ElmtID = Array<OneD, int>(nbcs,-1);
1406  }
1407  else
1408  {
1409  fill(ElmtID.get(), ElmtID.get()+nbcs, -1);
1410  }
1411 
1412  if (VertID.num_elements() != nbcs)
1413  {
1414  VertID = Array<OneD, int>(nbcs);
1415  }
1416 
1417  // setup map of all global ids along boundary
1418  for (cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
1419  {
1420  Vid = m_bndCondExpansions[n]->GetExp(0)->GetGeom()->GetVertex(0)->GetVid();
1421  VertGID[Vid] = cnt++;
1422  }
1423 
1424  // Loop over elements and find verts that match;
1426  for(cnt = n = 0; n < GetExpSize(); ++n)
1427  {
1428  exp = (*m_exp)[n];
1429  for(i = 0; i < exp->GetNverts(); ++i)
1430  {
1431  id = exp->GetGeom()->GetVid(i);
1432 
1433  if (VertGID.count(id) > 0)
1434  {
1435  bid = VertGID.find(id)->second;
1436  ASSERTL1(ElmtID[bid] == -1,"Edge already set");
1437  ElmtID[bid] = n;
1438  VertID[bid] = i;
1439  cnt ++;
1440  }
1441  }
1442  }
1443 
1444  ASSERTL1(cnt == nbcs,"Failed to visit all boundary condtiions");
1445  }
1446 
1448  boost::shared_ptr<ExpList> &result)
1449  {
1450  int n, cnt, nq;
1451  int offsetOld, offsetNew;
1452  Array<OneD, NekDouble> tmp1, tmp2;
1453  std::vector<unsigned int> eIDs;
1454 
1455  Array<OneD, int> ElmtID,EdgeID;
1456  GetBoundaryToElmtMap(ElmtID,EdgeID);
1457 
1458  // Skip other boundary regions
1459  for (cnt = n = 0; n < i; ++n)
1460  {
1461  cnt += m_bndCondExpansions[n]->GetExpSize();
1462  }
1463 
1464  // Populate eIDs with information from BoundaryToElmtMap
1465  for (n = 0; n < m_bndCondExpansions[i]->GetExpSize(); ++n)
1466  {
1467  eIDs.push_back(ElmtID[cnt+n]);
1468  }
1469 
1470  // Create expansion list
1471  result =
1473 
1474  // Copy phys and coeffs to new explist
1475  for (n = 0; n < result->GetExpSize(); ++n)
1476  {
1477  nq = GetExp(ElmtID[cnt+n])->GetTotPoints();
1478  offsetOld = GetPhys_Offset(ElmtID[cnt+n]);
1479  offsetNew = result->GetPhys_Offset(n);
1480  Vmath::Vcopy(nq, tmp1 = GetPhys()+ offsetOld, 1,
1481  tmp2 = result->UpdatePhys()+ offsetNew, 1);
1482 
1483  nq = GetExp(ElmtID[cnt+n])->GetNcoeffs();
1484  offsetOld = GetCoeff_Offset(ElmtID[cnt+n]);
1485  offsetNew = result->GetCoeff_Offset(n);
1486  Vmath::Vcopy(nq, tmp1 = GetCoeffs()+ offsetOld, 1,
1487  tmp2 = result->UpdateCoeffs()+ offsetNew, 1);
1488  }
1489  }
1490 
1491  /**
1492  * @brief Reset this field, so that geometry information can be updated.
1493  */
1495  {
1496  ExpList::v_Reset();
1497 
1498  // Reset boundary condition expansions.
1499  for (int n = 0; n < m_bndCondExpansions.num_elements(); ++n)
1500  {
1501  m_bndCondExpansions[n]->Reset();
1502  }
1503  }
1504 
1505  /**
1506  * Search through the edge expansions and identify which ones
1507  * have Robin/Mixed type boundary conditions. If find a Robin
1508  * boundary then store the edge id of the boundary condition
1509  * and the array of points of the physical space boundary
1510  * condition which are hold the boundary condition primitive
1511  * variable coefficient at the quatrature points
1512  *
1513  * \return std map containing the robin boundary condition
1514  * info using a key of the element id
1515  *
1516  * There is a next member to allow for more than one Robin
1517  * boundary condition per element
1518  */
1519  map<int, RobinBCInfoSharedPtr> DisContField1D::v_GetRobinBCInfo(void)
1520  {
1521  int i;
1522  map<int, RobinBCInfoSharedPtr> returnval;
1523  Array<OneD, int> ElmtID,VertID;
1524  GetBoundaryToElmtMap(ElmtID,VertID);
1525 
1526  for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
1527  {
1528  MultiRegions::ExpListSharedPtr locExpList;
1529 
1530  if (m_bndConditions[i]->GetBoundaryConditionType() ==
1532  {
1533  int elmtid;
1534  Array<OneD, NekDouble> Array_tmp;
1535 
1536  locExpList = m_bndCondExpansions[i];
1537 
1538  RobinBCInfoSharedPtr rInfo =
1541  VertID[i],Array_tmp = locExpList->GetPhys());
1542 
1543  elmtid = ElmtID[i];
1544  // make link list if necessary (not likely in
1545  // 1D but needed in 2D & 3D)
1546  if(returnval.count(elmtid) != 0)
1547  {
1548  rInfo->next = returnval.find(elmtid)->second;
1549  }
1550  returnval[elmtid] = rInfo;
1551  }
1552  }
1553 
1554  return returnval;
1555  }
1556  } // end of namespace
1557 } //end of namespace
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:814
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1837
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:1437
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n...
Definition: ExpList.h:1929
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
Evaluate all boundary conditions at a given time..
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1398
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2094
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
bool isLocal
Flag specifying if this entity is local to this partition.
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n...
Definition: ExpList.h:1937
const BoundaryConditionCollection & GetBoundaryConditions(void) const
Definition: Conditions.h:237
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
boost::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo()
void FindPeriodicVertices(const SpatialDomains::BoundaryConditions &bcs, const std::string variable)
Generate a associative map of periodic vertices in a mesh.
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1920
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
STL namespace.
SpatialDomains::BoundaryConditionsSharedPtr GetDomainBCs(const SpatialDomains::CompositeMap &domain, const SpatialDomains::BoundaryConditions &Allbcs, const std::string &variable)
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
void SetUpDG(const std::string &variable)
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2147
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:1424
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
boost::shared_ptr< Expansion0D > Expansion0DSharedPtr
Definition: Expansion0D.h:49
boost::shared_ptr< GlobalLinSys > GenGlobalBndLinSys(const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocTo...
Definition: ExpList.cpp:1194
void SetAdjacentElementExp(int vertex, Expansion1DSharedPtr &v)
Definition: Expansion0D.h:115
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:47
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
BasisManagerT & BasisManager(void)
GlobalLinSysSharedPtr GetGlobalBndLinSys(const GlobalLinSysKey &mkey)
For a given key, returns the associated global linear system.
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:988
void SetBoundaryConditionExpansion(const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::BoundaryConditions &bcs, const std::string variable, Array< OneD, MultiRegions::ExpListSharedPtr > &bndCondExpansions, Array< OneD, SpatialDomains::BoundaryConditionShPtr > &bndConditions)
Populates the list of boundary condition expansions.
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
The base class for all shapes.
Definition: StdExpansion.h:69
boost::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2004
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
virtual void v_GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result)
DisContField1D()
Default constructor.
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
int id
Geometry ID of entity.
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:965
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
std::map< int, Composite >::const_iterator CompositeMapConstIter
Definition: MeshGraph.h:117
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:206
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
This class is the abstraction of a global discontinuous two- dimensional spectral/hp element expansio...
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:913
void GenerateBoundaryConditionExpansion(const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::BoundaryConditions &bcs, const std::string variable)
Discretises the boundary conditions.
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
AssemblyMapDGSharedPtr m_traceMap
Local to global DG mapping for trace space.
Defines a specification for a set of points.
Definition: Points.h:58
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
double NekDouble
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
Defines a list of flags.
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:215
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
virtual void v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
Solve the Helmholtz equation.
static const NekDouble kNekUnsetDouble
Describe a linear system.
std::map< int, Composite > CompositeMap
Definition: MeshGraph.h:115
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
Describes a matrix with ordering defined by a local to global map.
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:1890
bool IsLeftAdjacentVertex(const int n, const int e)
GlobalLinSysMapShPtr m_globalBndMat
Global boundary matrix.
boost::shared_ptr< ExpList0D > ExpList0DSharedPtr
Shared pointer to an ExpList0D object.
Definition: ExpList0D.h:54
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
Definition: AssemblyMap.h:55
boost::shared_ptr< BoundaryConditionMap > BoundaryConditionMapShPtr
Definition: Conditions.h:214
std::vector< bool > & GetNegatedFluxNormal(void)
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
Discretised boundary conditions.
virtual void v_Reset()
Reset this field, so that geometry information can be updated.
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:1868
boost::shared_ptr< MeshGraph1D > MeshGraph1DSharedPtr
Definition: MeshGraph1D.h:86
std::set< int > m_boundaryVerts
A set storing the global IDs of any boundary edges.
boost::shared_ptr< BoundaryConditions > BoundaryConditionsSharedPtr
Definition: Conditions.h:271
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1562
This class is the abstraction of a one-dimensional multi-elemental expansions which is merely a colle...
Definition: ExpList1D.h:61
ExpListSharedPtr m_trace
Trace space storage for points between elements.
std::vector< int > m_periodicFwdCopy
A vector indicating degress of freedom which need to be copied from forwards to backwards space in ca...
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:2864
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:52
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:227
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:208
boost::shared_ptr< Basis > BasisSharedPtr
PeriodicMap m_periodicVerts
A map which identifies groups of periodic vertices.
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:723
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
boost::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:53
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
boost::shared_ptr< BoundaryRegion > BoundaryRegionShPtr
Definition: Conditions.h:204
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &VertID)