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