Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DisContField2D.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 //
3 // File DisContField2D.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 2D domain with boundary conditions using
33 // LDG flux.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 #include <LocalRegions/MatrixKey.h>
40 #include <LocalRegions/Expansion.h>
41 #include <LocalRegions/QuadExp.h>
42 #include <LocalRegions/TriExp.h>
46 
47 using namespace std;
48 
49 namespace Nektar
50 {
51  namespace MultiRegions
52  {
53  /**
54  * @class DisContField2D
55  * Abstraction of a global discontinuous two-dimensional spectral/hp
56  * element expansion which approximates the solution of a set of
57  * partial differential equations.
58  */
59 
60  /**
61  * @brief Default constructor.
62  */
63  DisContField2D::DisContField2D(void)
64  : ExpList2D (),
65  m_bndCondExpansions(),
66  m_bndConditions (),
67  m_trace (NullExpListSharedPtr)
68  {
69  }
70 
72  const DisContField2D &In,
73  const bool DeclareCoeffPhysArrays)
74  : ExpList2D (In,DeclareCoeffPhysArrays),
75  m_bndCondExpansions (In.m_bndCondExpansions),
76  m_bndConditions (In.m_bndConditions),
77  m_globalBndMat (In.m_globalBndMat),
78  m_traceMap (In.m_traceMap),
79  m_locTraceToTraceMap (In.m_locTraceToTraceMap),
80  m_boundaryEdges (In.m_boundaryEdges),
81  m_periodicVerts (In.m_periodicVerts),
82  m_periodicEdges (In.m_periodicEdges),
83  m_periodicFwdCopy (In.m_periodicFwdCopy),
84  m_periodicBwdCopy (In.m_periodicBwdCopy),
85  m_leftAdjacentEdges (In.m_leftAdjacentEdges)
86  {
87  if (In.m_trace)
88  {
90  *boost::dynamic_pointer_cast<ExpList1D>(In.m_trace),
91  DeclareCoeffPhysArrays);
92  }
93  }
94 
95  /**
96  * @brief Constructs a global discontinuous field based on an input
97  * mesh with boundary conditions.
98  */
100  const LibUtilities::SessionReaderSharedPtr &pSession,
101  const SpatialDomains::MeshGraphSharedPtr &graph2D,
102  const std::string &variable,
103  const bool SetUpJustDG,
104  const bool DeclareCoeffPhysArrays)
105  : ExpList2D(pSession, graph2D, DeclareCoeffPhysArrays, variable),
106  m_bndCondExpansions(),
107  m_bndConditions(),
108  m_trace(NullExpListSharedPtr),
109  m_periodicVerts(),
110  m_periodicEdges(),
111  m_periodicFwdCopy(),
112  m_periodicBwdCopy()
113  {
114 
115  if (variable.compare("DefaultVar") != 0) // do not set up BCs if default variable
116  {
118  GenerateBoundaryConditionExpansion(graph2D, bcs, variable,
119  DeclareCoeffPhysArrays);
120 
121  if (DeclareCoeffPhysArrays)
122  {
123  EvaluateBoundaryConditions(0.0, variable);
124  }
125 
126  // Find periodic edges for this variable.
127  FindPeriodicEdges(bcs, variable);
128  }
129 
130  if (SetUpJustDG)
131  {
132  SetUpDG(variable);
133  }
134  else
135  {
136  // Set element edges to point to Robin BC edges if required.
137  int i, cnt;
138  Array<OneD, int> ElmtID, EdgeID;
139  GetBoundaryToElmtMap(ElmtID, EdgeID);
140 
141  for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
142  {
144  int e;
145  locExpList = m_bndCondExpansions[i];
146 
147  for(e = 0; e < locExpList->GetExpSize(); ++e)
148  {
150  (*m_exp)[ElmtID[cnt+e]]->
151  as<LocalRegions::Expansion2D>();
153  locExpList->GetExp(e)->
154  as<LocalRegions::Expansion1D>();
156  locExpList->GetExp(e)->
157  as<LocalRegions::Expansion> ();
158 
159  exp2d->SetEdgeExp(EdgeID[cnt+e], exp1d);
160  exp1d->SetAdjacentElementExp(EdgeID[cnt+e], exp2d);
161  }
162  cnt += m_bndCondExpansions[i]->GetExpSize();
163  }
164 
165  if(m_session->DefinesSolverInfo("PROJECTION"))
166  {
167  std::string ProjectStr =
168  m_session->GetSolverInfo("PROJECTION");
169  if((ProjectStr == "MixedCGDG") ||
170  (ProjectStr == "Mixed_CG_Discontinuous"))
171  {
172  SetUpDG();
173  }
174  else
175  {
177  }
178  }
179  else
180  {
182  }
183  }
184  }
185 
186 
187  /*
188  * @brief Copy type constructor which declares new boundary conditions
189  * and re-uses mapping info and trace space if possible
190  */
192  const DisContField2D &In,
193  const SpatialDomains::MeshGraphSharedPtr &graph2D,
194  const std::string &variable,
195  const bool SetUpJustDG,
196  const bool DeclareCoeffPhysArrays)
197  : ExpList2D(In,DeclareCoeffPhysArrays),
198  m_trace(NullExpListSharedPtr)
199  {
200  // Set up boundary conditions for this variable.
201  // Do not set up BCs if default variable
202  if(variable.compare("DefaultVar") != 0)
203  {
205  GenerateBoundaryConditionExpansion(graph2D, bcs, variable);
206 
207  if (DeclareCoeffPhysArrays)
208  {
209  EvaluateBoundaryConditions(0.0, variable);
210  }
211 
213  {
214  // Find periodic edges for this variable.
215  FindPeriodicEdges(bcs, variable);
216 
217  if(SetUpJustDG)
218  {
219  SetUpDG();
220  }
221  else
222  {
223  // set elmt edges to point to robin bc edges if required
224  int i, cnt = 0;
225  Array<OneD, int> ElmtID,EdgeID;
226  GetBoundaryToElmtMap(ElmtID,EdgeID);
227 
228  for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
229  {
231 
232  int e;
233  locExpList = m_bndCondExpansions[i];
234 
235  for(e = 0; e < locExpList->GetExpSize(); ++e)
236  {
238  = (*m_exp)[ElmtID[cnt+e]]->
239  as<LocalRegions::Expansion2D>();
241  = locExpList->GetExp(e)->
242  as<LocalRegions::Expansion1D>();
244  = locExpList->GetExp(e)->
245  as<LocalRegions::Expansion> ();
246 
247  exp2d->SetEdgeExp(EdgeID[cnt+e], exp1d);
248  exp1d->SetAdjacentElementExp(EdgeID[cnt+e],
249  exp2d);
250  }
251 
252  cnt += m_bndCondExpansions[i]->GetExpSize();
253  }
254 
255 
256  if (m_session->DefinesSolverInfo("PROJECTION"))
257  {
258  std::string ProjectStr =
259  m_session->GetSolverInfo("PROJECTION");
260 
261  if ((ProjectStr == "MixedCGDG") ||
262  (ProjectStr == "Mixed_CG_Discontinuous"))
263  {
264  SetUpDG();
265  }
266  else
267  {
269  }
270  }
271  else
272  {
274  }
275  }
276  }
277  else
278  {
279  if (SetUpJustDG)
280  {
282  m_trace = In.m_trace;
283  m_traceMap = In.m_traceMap;
291  }
292  else
293  {
295  m_trace = In.m_trace;
296  m_traceMap = In.m_traceMap;
304 
305  // set elmt edges to point to robin bc edges if required
306  int i, cnt = 0;
307  Array<OneD, int> ElmtID, EdgeID;
308  GetBoundaryToElmtMap(ElmtID, EdgeID);
309 
310  for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
311  {
313 
314  int e;
315  locExpList = m_bndCondExpansions[i];
316 
317  for (e = 0; e < locExpList->GetExpSize(); ++e)
318  {
320  = (*m_exp)[ElmtID[cnt+e]]->
321  as<LocalRegions::Expansion2D>();
323  = locExpList->GetExp(e)->
324  as<LocalRegions::Expansion1D>();
326  = locExpList->GetExp(e)->
327  as<LocalRegions::Expansion> ();
328 
329  exp2d->SetEdgeExp(EdgeID[cnt+e], exp1d);
330  exp1d->SetAdjacentElementExp(EdgeID[cnt+e],
331  exp2d);
332  }
333  cnt += m_bndCondExpansions[i]->GetExpSize();
334  }
335 
337  }
338  }
339  }
340  }
341 
342  /**
343  * @brief Default destructor.
344  */
346  {
347  }
348 
350  const GlobalLinSysKey &mkey)
351  {
353  "Routine currently only tested for HybridDGHelmholtz");
354  ASSERTL1(mkey.GetGlobalSysSolnType() ==
355  m_traceMap->GetGlobalSysSolnType(),
356  "The local to global map is not set up for the requested "
357  "solution type");
358 
359  GlobalLinSysSharedPtr glo_matrix;
360  GlobalLinSysMap::iterator matrixIter = m_globalBndMat->find(mkey);
361 
362  if(matrixIter == m_globalBndMat->end())
363  {
364  glo_matrix = GenGlobalBndLinSys(mkey,m_traceMap);
365  (*m_globalBndMat)[mkey] = glo_matrix;
366  }
367  else
368  {
369  glo_matrix = matrixIter->second;
370  }
371 
372  return glo_matrix;
373  }
374 
375  /**
376  * @brief Set up all DG member variables and maps.
377  */
378  void DisContField2D::SetUpDG(const std::string variable)
379  {
380  // Check for multiple calls
382  {
383  return;
384  }
385 
386  ExpList1DSharedPtr trace;
388  boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(
389  m_graph);
390 
391  // Set up matrix map
394 
395  // Set up trace space
398  graph2D, m_periodicEdges);
399 
400  m_trace = boost::dynamic_pointer_cast<ExpList>(trace);
402  AllocateSharedPtr(m_session, graph2D, trace, *this,
405  variable);
406 
407  if (m_session->DefinesCmdLineArgument("verbose"))
408  {
409  m_traceMap->PrintStats(std::cout, variable);
410  }
411 
413  &elmtToTrace = m_traceMap->GetElmtToTrace();
414 
415  // Scatter trace segments to 2D elements. For each element,
416  // we find the trace segment associated to each edge. The
417  // element then retains a pointer to the trace space segments,
418  // to ensure uniqueness of normals when retrieving from two
419  //adjoining elements which do not lie in a plane.
420  for (int i = 0; i < m_exp->size(); ++i)
421  {
422  for (int j = 0; j < (*m_exp)[i]->GetNedges(); ++j)
423  {
425  (*m_exp)[i]->as<LocalRegions::Expansion2D>();
427  elmtToTrace[i][j]->as<LocalRegions::Expansion1D>();
428  LocalRegions::ExpansionSharedPtr exp = elmtToTrace[i][j];;
429  exp2d->SetEdgeExp (j, exp1d);
430  exp1d->SetAdjacentElementExp(j, exp2d);
431  }
432  }
433 
434  // Set up physical normals
436 
437  // Set up information for parallel and periodic problems.
438  for (int i = 0; i < m_trace->GetExpSize(); ++i)
439  {
441  m_trace->GetExp(i)->as<LocalRegions::Expansion1D>();
442 
443  int offset = m_trace->GetPhys_Offset(i);
444  int traceGeomId = traceEl->GetGeom1D()->GetGlobalID();
445  PeriodicMap::iterator pIt = m_periodicEdges.find(traceGeomId);
446 
447  if (pIt != m_periodicEdges.end() && !pIt->second[0].isLocal)
448  {
449  if (traceGeomId != min(pIt->second[0].id, traceGeomId))
450  {
451  traceEl->GetLeftAdjacentElementExp()->NegateEdgeNormal(
452  traceEl->GetLeftAdjacentElementEdge());
453  }
454  }
455  else if (m_traceMap->GetTraceToUniversalMapUnique(offset) < 0)
456  {
457  traceEl->GetLeftAdjacentElementExp()->NegateEdgeNormal(
458  traceEl->GetLeftAdjacentElementEdge());
459  }
460  }
461 
462  int cnt, n, e;
463 
464  // Identify boundary edges
465  for (cnt = 0, n = 0; n < m_bndCondExpansions.num_elements(); ++n)
466  {
467  if (m_bndConditions[n]->GetBoundaryConditionType() !=
469  {
470  for (e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
471  {
472  m_boundaryEdges.insert(
473  m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
474  }
475  }
476  cnt += m_bndCondExpansions[n]->GetExpSize();
477  }
478 
479  // Set up information for periodic boundary conditions.
480  boost::unordered_map<int,pair<int,int> > perEdgeToExpMap;
481  boost::unordered_map<int,pair<int,int> >::iterator it2;
482  for (cnt = n = 0; n < m_exp->size(); ++n)
483  {
484  for (e = 0; e < (*m_exp)[n]->GetNedges(); ++e, ++cnt)
485  {
487  (*m_exp)[n]->GetGeom()->GetEid(e));
488 
489  if (it != m_periodicEdges.end())
490  {
491  perEdgeToExpMap[it->first] = make_pair(n, e);
492  }
493  }
494  }
495 
496  // Set up left-adjacent edge list.
497  m_leftAdjacentEdges.resize(cnt);
498  cnt = 0;
499  for (int i = 0; i < m_exp->size(); ++i)
500  {
501  for (int j = 0; j < (*m_exp)[i]->GetNedges(); ++j, ++cnt)
502  {
504  }
505  }
506 
507  // Set up mapping to copy Fwd of periodic bcs to Bwd of other edge.
508  cnt = 0;
509  for (int n = 0; n < m_exp->size(); ++n)
510  {
511  for (int e = 0; e < (*m_exp)[n]->GetNedges(); ++e, ++cnt)
512  {
513  int edgeGeomId = (*m_exp)[n]->GetGeom()->GetEid(e);
514  int offset = m_trace->GetPhys_Offset(
515  elmtToTrace[n][e]->GetElmtId());
516 
517  // Check to see if this face is periodic.
518  PeriodicMap::iterator it = m_periodicEdges.find(edgeGeomId);
519 
520  if (it != m_periodicEdges.end())
521  {
522  const PeriodicEntity &ent = it->second[0];
523  it2 = perEdgeToExpMap.find(ent.id);
524 
525  if (it2 == perEdgeToExpMap.end())
526  {
527  if (m_session->GetComm()->
528  GetRowComm()->GetSize() > 1 && !ent.isLocal)
529  {
530  continue;
531  }
532  else
533  {
534  ASSERTL1(false, "Periodic edge not found!");
535  }
536  }
537 
539  "Periodic edge in non-forward space?");
540 
541  int offset2 = m_trace->GetPhys_Offset(
542  elmtToTrace[it2->second.first][it2->second.second]->
543  GetElmtId());
544 
545  // Calculate relative orientations between edges to
546  // calculate copying map.
547  int nquad = elmtToTrace[n][e]->GetNumPoints(0);
548 
549  vector<int> tmpBwd(nquad);
550  vector<int> tmpFwd(nquad);
551 
552  if (ent.orient == StdRegions::eForwards)
553  {
554  for (int i = 0; i < nquad; ++i)
555  {
556  tmpBwd[i] = offset2 + i;
557  tmpFwd[i] = offset + i;
558  }
559  }
560  else
561  {
562  for (int i = 0; i < nquad; ++i)
563  {
564  tmpBwd[i] = offset2 + i;
565  tmpFwd[i] = offset + nquad - i - 1;
566  }
567  }
568 
569  for (int i = 0; i < nquad; ++i)
570  {
571  m_periodicFwdCopy.push_back(tmpFwd[i]);
572  m_periodicBwdCopy.push_back(tmpBwd[i]);
573  }
574  }
575  }
576  }
577 
579  AllocateSharedPtr(*this, m_trace, elmtToTrace,
581  }
582 
583  /**
584  * For each boundary region, checks that the types and number of
585  * boundary expansions in that region match.
586  * @param In ContField2D to compare with.
587  * @return True if boundary conditions match.
588  */
590  const DisContField2D &In)
591  {
592  int i;
593  bool returnval = true;
594 
595  for(i = 0; i < m_bndConditions.num_elements(); ++i)
596  {
597 
598  // check to see if boundary condition type is the same
599  // and there are the same number of boundary
600  // conditions in the boundary definition.
601  if((m_bndConditions[i]->GetBoundaryConditionType()
602  != In.m_bndConditions[i]->GetBoundaryConditionType())||
604  != In.m_bndCondExpansions[i]->GetExpSize()))
605  {
606  returnval = false;
607  break;
608  }
609  }
610 
611  // Compare with all other processes. Return true only if all
612  // processes report having the same boundary conditions.
613  int vSame = (returnval?1:0);
614  m_comm->GetRowComm()->AllReduce(vSame, LibUtilities::ReduceMin);
615 
616  return (vSame == 1);
617  }
618 
619  /**
620  * \brief This function discretises the boundary conditions by setting
621  * up a list of one-dimensional boundary expansions.
622  *
623  * According to their boundary region, the separate segmental boundary
624  * expansions are bundled together in an object of the class
625  * MultiRegions#ExpList1D.
626  *
627  * \param graph2D A mesh, containing information about the domain and
628  * the spectral/hp element expansion.
629  * \param bcs An entity containing information about the
630  * boundaries and boundary conditions.
631  * \param variable An optional parameter to indicate for which variable
632  * the boundary conditions should be discretised.
633  */
635  const SpatialDomains::MeshGraphSharedPtr &graph2D,
637  const std::string &variable,
638  const bool DeclareCoeffPhysArrays)
639  {
640  int cnt = 0;
644  bcs.GetBoundaryRegions();
645  const SpatialDomains::BoundaryConditionCollection &bconditions =
646  bcs.GetBoundaryConditions();
647  SpatialDomains::BoundaryRegionCollection::const_iterator it;
648 
649  // count the number of non-periodic boundary regions
650  for (it = bregions.begin(); it != bregions.end(); ++it)
651  {
652  bc = GetBoundaryCondition(bconditions, it->first, variable);
653 
654  if (bc->GetBoundaryConditionType() != SpatialDomains::ePeriodic)
655  {
656  cnt++;
657  }
658  }
659 
662  m_bndConditions =
664 
665  cnt = 0;
666 
667  // list non-periodic boundaries
668  for (it = bregions.begin(); it != bregions.end(); ++it)
669  {
670  bc = GetBoundaryCondition(bconditions, it->first, variable);
671 
672  if (bc->GetBoundaryConditionType() != SpatialDomains::ePeriodic)
673  {
675  ::AllocateSharedPtr(m_session, *(it->second), graph2D,
676  DeclareCoeffPhysArrays, variable);
677 
678  m_bndCondExpansions[cnt] = locExpList;
679  m_bndConditions[cnt] = bc;
680 
681 
682  std::string type = m_bndConditions[cnt]->GetUserDefined();
683 
684  // Set up normals on non-Dirichlet boundary
685  // conditions. Second two conditions ideally
686  // should be in local solver setup (when made into factory)
687  if((bc->GetBoundaryConditionType() !=
689  boost::iequals(type,"I") ||
690  boost::iequals(type,"CalcBC"))
691  {
693  }
694 
695  cnt++;
696  }
697  }
698  }
699 
700  /**
701  * @brief Determine the periodic edges and vertices for the given graph.
702  *
703  * Note that much of this routine is the same as the three-dimensional
704  * version, which therefore has much better documentation.
705  *
706  * @param bcs Information about the boundary conditions.
707  * @param variable Specifies the field.
708  *
709  * @see DisContField3D::FindPeriodicFaces
710  */
713  const std::string &variable)
714  {
716  = bcs.GetBoundaryRegions();
718  = bcs.GetBoundaryConditions();
720  = boost::dynamic_pointer_cast<
722  SpatialDomains::BoundaryRegionCollection::const_iterator it;
723 
725  m_session->GetComm()->GetRowComm();
727  m_session->GetCompositeOrdering();
728  LibUtilities::BndRegionOrdering bndRegOrder =
729  m_session->GetBndRegionOrdering();
731  m_graph->GetComposites();
732 
733  // Unique collection of pairs of periodic composites (i.e. if
734  // composites 1 and 2 are periodic then this map will contain either
735  // the pair (1,2) or (2,1) but not both).
736  map<int,int> perComps;
737  map<int,vector<int> > allVerts;
738  set<int> locVerts;
739  map<int,StdRegions::Orientation> allEdges;
740 
741  int region1ID, region2ID, i, j, k, cnt;
743 
744  // Set up a set of all local verts and edges.
745  for(i = 0; i < (*m_exp).size(); ++i)
746  {
747  for(j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
748  {
749  int id = (*m_exp)[i]->GetGeom()->GetVid(j);
750  locVerts.insert(id);
751  }
752  }
753 
754  // Construct list of all periodic pairs local to this process.
755  for (it = bregions.begin(); it != bregions.end(); ++it)
756  {
757  locBCond = GetBoundaryCondition(
758  bconditions, it->first, variable);
759 
760  if (locBCond->GetBoundaryConditionType()
762  {
763  continue;
764  }
765 
766  // Identify periodic boundary region IDs.
767  region1ID = it->first;
768  region2ID = boost::static_pointer_cast<
770  locBCond)->m_connectedBoundaryRegion;
771 
772  // From this identify composites. Note that in serial this will
773  // be an empty map.
774  int cId1, cId2;
775  if (vComm->GetSize() == 1)
776  {
777  cId1 = it->second->begin()->first;
778  cId2 = bregions.find(region2ID)->second->begin()->first;
779  }
780  else
781  {
782  cId1 = bndRegOrder.find(region1ID)->second[0];
783  cId2 = bndRegOrder.find(region2ID)->second[0];
784  }
785 
786  ASSERTL0(it->second->size() == 1,
787  "Boundary region "+boost::lexical_cast<string>(
788  region1ID)+" should only contain 1 composite.");
789 
790  // Construct set containing all periodic edges on this process
791  SpatialDomains::Composite c = it->second->begin()->second;
792 
793  vector<unsigned int> tmpOrder;
794 
795  for (i = 0; i < c->size(); ++i)
796  {
798  boost::dynamic_pointer_cast<
799  SpatialDomains::SegGeom>((*c)[i]);
800  ASSERTL0(segGeom, "Unable to cast to shared ptr");
801 
803  graph2D->GetElementsFromEdge(segGeom);
804  ASSERTL0(elmt->size() == 1,
805  "The periodic boundaries belong to "
806  "more than one element of the mesh");
807 
809  boost::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
810  (*elmt)[0]->m_Element);
811 
812  allEdges[(*c)[i]->GetGlobalID()] =
813  geom->GetEorient((*elmt)[0]->m_EdgeIndx);
814 
815  // In serial mesh partitioning will not have occurred so
816  // need to fill composite ordering map manually.
817  if (vComm->GetSize() == 1)
818  {
819  tmpOrder.push_back((*c)[i]->GetGlobalID());
820  }
821 
822  vector<int> vertList(2);
823  vertList[0] = segGeom->GetVid(0);
824  vertList[1] = segGeom->GetVid(1);
825  allVerts[(*c)[i]->GetGlobalID()] = vertList;
826  }
827 
828  if (vComm->GetSize() == 1)
829  {
830  compOrder[it->second->begin()->first] = tmpOrder;
831  }
832 
833  // See if we already have either region1 or region2 stored in
834  // perComps map.
835  if (perComps.count(cId1) == 0)
836  {
837  if (perComps.count(cId2) == 0)
838  {
839  perComps[cId1] = cId2;
840  }
841  else
842  {
843  std::stringstream ss;
844  ss << "Boundary region " << cId2 << " should be "
845  << "periodic with " << perComps[cId2] << " but "
846  << "found " << cId1 << " instead!";
847  ASSERTL0(perComps[cId2] == cId1, ss.str());
848  }
849  }
850  else
851  {
852  std::stringstream ss;
853  ss << "Boundary region " << cId1 << " should be "
854  << "periodic with " << perComps[cId1] << " but "
855  << "found " << cId2 << " instead!";
856  ASSERTL0(perComps[cId1] == cId1, ss.str());
857  }
858  }
859 
860  // Process local edge list to obtain relative edge orientations.
861  int n = vComm->GetSize();
862  int p = vComm->GetRank();
863  int totEdges;
864  Array<OneD, int> edgecounts(n,0);
865  Array<OneD, int> edgeoffset(n,0);
866  Array<OneD, int> vertoffset(n,0);
867 
868  edgecounts[p] = allEdges.size();
869  vComm->AllReduce(edgecounts, LibUtilities::ReduceSum);
870 
871  edgeoffset[0] = 0;
872  for (i = 1; i < n; ++i)
873  {
874  edgeoffset[i] = edgeoffset[i-1] + edgecounts[i-1];
875  }
876 
877  totEdges = Vmath::Vsum(n, edgecounts, 1);
878  Array<OneD, int> edgeIds (totEdges, 0);
879  Array<OneD, int> edgeOrient(totEdges, 0);
880  Array<OneD, int> edgeVerts (totEdges, 0);
881 
883 
884  for (i = 0, sIt = allEdges.begin(); sIt != allEdges.end(); ++sIt)
885  {
886  edgeIds [edgeoffset[p] + i ] = sIt->first;
887  edgeOrient[edgeoffset[p] + i ] = sIt->second;
888  edgeVerts [edgeoffset[p] + i++] = allVerts[sIt->first].size();
889  }
890 
891  vComm->AllReduce(edgeIds, LibUtilities::ReduceSum);
892  vComm->AllReduce(edgeOrient, LibUtilities::ReduceSum);
893  vComm->AllReduce(edgeVerts, LibUtilities::ReduceSum);
894 
895  // Calculate number of vertices on each processor.
896  Array<OneD, int> procVerts(n,0);
897  int nTotVerts;
898 
899  // Note if there are no periodic edges at all calling Vsum will
900  // cause a segfault.
901  if (totEdges > 0)
902  {
903  nTotVerts = Vmath::Vsum(totEdges, edgeVerts, 1);
904  }
905  else
906  {
907  nTotVerts = 0;
908  }
909 
910  for (i = 0; i < n; ++i)
911  {
912  if (edgecounts[i] > 0)
913  {
914  procVerts[i] = Vmath::Vsum(
915  edgecounts[i], edgeVerts + edgeoffset[i], 1);
916  }
917  else
918  {
919  procVerts[i] = 0;
920  }
921  }
922  vertoffset[0] = 0;
923 
924  for (i = 1; i < n; ++i)
925  {
926  vertoffset[i] = vertoffset[i-1] + procVerts[i-1];
927  }
928 
929  Array<OneD, int> vertIds(nTotVerts, 0);
930  for (i = 0, sIt = allEdges.begin(); sIt != allEdges.end(); ++sIt)
931  {
932  for (j = 0; j < allVerts[sIt->first].size(); ++j)
933  {
934  vertIds[vertoffset[p] + i++] = allVerts[sIt->first][j];
935  }
936  }
937 
938  vComm->AllReduce(vertIds, LibUtilities::ReduceSum);
939 
940  // For simplicity's sake create a map of edge id -> orientation.
941  map<int, StdRegions::Orientation> orientMap;
942  map<int, vector<int> > vertMap;
943 
944  for (cnt = i = 0; i < totEdges; ++i)
945  {
946  ASSERTL0(orientMap.count(edgeIds[i]) == 0,
947  "Already found edge in orientation map!");
948  orientMap[edgeIds[i]] = (StdRegions::Orientation)edgeOrient[i];
949 
950  vector<int> verts(edgeVerts[i]);
951 
952  for (j = 0; j < edgeVerts[i]; ++j)
953  {
954  verts[j] = vertIds[cnt++];
955  }
956  vertMap[edgeIds[i]] = verts;
957  }
958 
959  // Go through list of composites and figure out which edges are
960  // parallel from original ordering in session file. This includes
961  // composites which are not necessarily on this process.
962  map<int,int>::iterator cIt, pIt;
963  map<int,int>::const_iterator oIt;
964 
965  map<int,int> allCompPairs;
966 
967  // Store temporary map of periodic vertices which will hold all
968  // periodic vertices on the entire mesh so that doubly periodic
969  // vertices can be counted properly across partitions. Local
970  // vertices are copied into m_periodicVerts at the end of the
971  // function.
972  PeriodicMap periodicVerts;
973 
974  for (cIt = perComps.begin(); cIt != perComps.end(); ++cIt)
975  {
977  const int id1 = cIt->first;
978  const int id2 = cIt->second;
979  std::string id1s = boost::lexical_cast<string>(id1);
980  std::string id2s = boost::lexical_cast<string>(id2);
981 
982  if (compMap.count(id1) > 0)
983  {
984  c[0] = compMap[id1];
985  }
986 
987  if (compMap.count(id2) > 0)
988  {
989  c[1] = compMap[id2];
990  }
991 
992  ASSERTL0(c[0] || c[1],
993  "Both composites not found on this process!");
994 
995  // Loop over composite ordering to construct list of all
996  // periodic edges regardless of whether they are on this
997  // process.
998  map<int,int> compPairs;
999 
1000  ASSERTL0(compOrder.count(id1) > 0,
1001  "Unable to find composite "+id1s+" in order map.");
1002  ASSERTL0(compOrder.count(id2) > 0,
1003  "Unable to find composite "+id2s+" in order map.");
1004  ASSERTL0(compOrder[id1].size() == compOrder[id2].size(),
1005  "Periodic composites "+id1s+" and "+id2s+
1006  " should have the same number of elements.");
1007  ASSERTL0(compOrder[id1].size() > 0,
1008  "Periodic composites "+id1s+" and "+id2s+
1009  " are empty!");
1010 
1011  // TODO: Add more checks.
1012  for (i = 0; i < compOrder[id1].size(); ++i)
1013  {
1014  int eId1 = compOrder[id1][i];
1015  int eId2 = compOrder[id2][i];
1016 
1017  ASSERTL0(compPairs.count(eId1) == 0,
1018  "Already paired.");
1019 
1020  if (compPairs.count(eId2) != 0)
1021  {
1022  ASSERTL0(compPairs[eId2] == eId1, "Pairing incorrect");
1023  }
1024  compPairs[eId1] = eId2;
1025  }
1026 
1027  // Construct set of all edges that we have locally on this
1028  // processor.
1029  set<int> locEdges;
1030  set<int>::iterator sIt;
1031  for (i = 0; i < 2; ++i)
1032  {
1033  if (!c[i])
1034  {
1035  continue;
1036  }
1037 
1038  if (c[i]->size() > 0)
1039  {
1040  for (j = 0; j < c[i]->size(); ++j)
1041  {
1042  locEdges.insert(c[i]->at(j)->GetGlobalID());
1043  }
1044  }
1045  }
1046 
1047  // Loop over all edges in the geometry composite.
1048  for (pIt = compPairs.begin(); pIt != compPairs.end(); ++pIt)
1049  {
1050  int ids [2] = {pIt->first, pIt->second};
1051  bool local[2] = {locEdges.count(pIt->first) > 0,
1052  locEdges.count(pIt->second) > 0};
1053 
1054  ASSERTL0(orientMap.count(ids[0]) > 0 &&
1055  orientMap.count(ids[1]) > 0,
1056  "Unable to find edge in orientation map.");
1057 
1058  allCompPairs[pIt->first ] = pIt->second;
1059  allCompPairs[pIt->second] = pIt->first;
1060 
1061  for (i = 0; i < 2; ++i)
1062  {
1063  if (!local[i])
1064  {
1065  continue;
1066  }
1067 
1068  int other = (i+1) % 2;
1069 
1071  orientMap[ids[i]] == orientMap[ids[other]] ?
1074 
1075  PeriodicEntity ent(ids [other], o,
1076  local[other]);
1077  m_periodicEdges[ids[i]].push_back(ent);
1078  }
1079 
1080  for (i = 0; i < 2; ++i)
1081  {
1082  int other = (i+1) % 2;
1083 
1085  orientMap[ids[i]] == orientMap[ids[other]] ?
1088 
1089  // Determine periodic vertices.
1090  vector<int> perVerts1 = vertMap[ids[i]];
1091  vector<int> perVerts2 = vertMap[ids[other]];
1092 
1093  map<int, pair<int, bool> > tmpMap;
1094  map<int, pair<int, bool> >::iterator mIt;
1095  if (o == StdRegions::eForwards)
1096  {
1097  tmpMap[perVerts1[0]] = make_pair(
1098  perVerts2[0], locVerts.count(perVerts2[0]) > 0);
1099  tmpMap[perVerts1[1]] = make_pair(
1100  perVerts2[1], locVerts.count(perVerts2[1]) > 0);
1101  }
1102  else
1103  {
1104  tmpMap[perVerts1[0]] = make_pair(
1105  perVerts2[1], locVerts.count(perVerts2[1]) > 0);
1106  tmpMap[perVerts1[1]] = make_pair(
1107  perVerts2[0], locVerts.count(perVerts2[0]) > 0);
1108  }
1109 
1110  for (mIt = tmpMap.begin(); mIt != tmpMap.end(); ++mIt)
1111  {
1112  // See if this vertex has been recorded already.
1113  PeriodicEntity ent2(mIt->second.first,
1115  mIt->second.second);
1116  PeriodicMap::iterator perIt = periodicVerts.find(
1117  mIt->first);
1118 
1119  if (perIt == periodicVerts.end())
1120  {
1121  periodicVerts[mIt->first].push_back(ent2);
1122  perIt = periodicVerts.find(mIt->first);
1123  }
1124  else
1125  {
1126  bool doAdd = true;
1127  for (j = 0; j < perIt->second.size(); ++j)
1128  {
1129  if (perIt->second[j].id == mIt->second.first)
1130  {
1131  doAdd = false;
1132  break;
1133  }
1134  }
1135 
1136  if (doAdd)
1137  {
1138  perIt->second.push_back(ent2);
1139  }
1140  }
1141  }
1142  }
1143  }
1144  }
1145 
1146  Array<OneD, int> pairSizes(n, 0);
1147  pairSizes[p] = allCompPairs.size();
1148  vComm->AllReduce(pairSizes, LibUtilities::ReduceSum);
1149 
1150  int totPairSizes = Vmath::Vsum(n, pairSizes, 1);
1151 
1152  Array<OneD, int> pairOffsets(n, 0);
1153  pairOffsets[0] = 0;
1154 
1155  for (i = 1; i < n; ++i)
1156  {
1157  pairOffsets[i] = pairOffsets[i-1] + pairSizes[i-1];
1158  }
1159 
1160  Array<OneD, int> first (totPairSizes, 0);
1161  Array<OneD, int> second(totPairSizes, 0);
1162 
1163  cnt = pairOffsets[p];
1164 
1165  for (pIt = allCompPairs.begin(); pIt != allCompPairs.end(); ++pIt)
1166  {
1167  first [cnt ] = pIt->first;
1168  second[cnt++] = pIt->second;
1169  }
1170 
1171  vComm->AllReduce(first, LibUtilities::ReduceSum);
1172  vComm->AllReduce(second, LibUtilities::ReduceSum);
1173 
1174  allCompPairs.clear();
1175 
1176  for(cnt = 0; cnt < totPairSizes; ++cnt)
1177  {
1178  allCompPairs[first[cnt]] = second[cnt];
1179  }
1180 
1181  // Search for periodic vertices and edges which are not in
1182  // a periodic composite but lie in this process. First, loop
1183  // over all information we have from other processors.
1184  for (cnt = i = 0; i < totEdges; ++i)
1185  {
1186  int edgeId = edgeIds[i];
1187 
1188  ASSERTL0(allCompPairs.count(edgeId) > 0,
1189  "Unable to find matching periodic edge.");
1190 
1191  int perEdgeId = allCompPairs[edgeId];
1192 
1193  for (j = 0; j < edgeVerts[i]; ++j, ++cnt)
1194  {
1195  int vId = vertIds[cnt];
1196 
1197  PeriodicMap::iterator perId = periodicVerts.find(vId);
1198 
1199  if (perId == periodicVerts.end())
1200  {
1201  // This vertex is not included in the map. Figure
1202  // out which vertex it is supposed to be periodic
1203  // with. perEdgeId is the edge ID which is periodic with
1204  // edgeId. The logic is much the same as the loop above.
1205  int perVertexId =
1206  orientMap[edgeId] == orientMap[perEdgeId] ?
1207  vertMap[perEdgeId][(j+1)%2] : vertMap[perEdgeId][j];
1208 
1209  PeriodicEntity ent(perVertexId,
1211  locVerts.count(perVertexId) > 0);
1212 
1213  periodicVerts[vId].push_back(ent);
1214  }
1215  }
1216  }
1217 
1218  // Loop over all periodic vertices to complete connectivity
1219  // information.
1220  PeriodicMap::iterator perIt, perIt2;
1221  for (perIt = periodicVerts.begin();
1222  perIt != periodicVerts.end(); ++perIt)
1223  {
1224  // Loop over associated vertices.
1225  for (i = 0; i < perIt->second.size(); ++i)
1226  {
1227  perIt2 = periodicVerts.find(perIt->second[i].id);
1228  ASSERTL0(perIt2 != periodicVerts.end(),
1229  "Couldn't find periodic vertex.");
1230 
1231  for (j = 0; j < perIt2->second.size(); ++j)
1232  {
1233  if (perIt2->second[j].id == perIt->first)
1234  {
1235  continue;
1236  }
1237 
1238  bool doAdd = true;
1239  for (k = 0; k < perIt->second.size(); ++k)
1240  {
1241  if (perIt2->second[j].id == perIt->second[k].id)
1242  {
1243  doAdd = false;
1244  break;
1245  }
1246  }
1247 
1248  if (doAdd)
1249  {
1250  perIt->second.push_back(perIt2->second[j]);
1251  }
1252  }
1253  }
1254  }
1255 
1256  // Do one final loop over periodic vertices to remove non-local
1257  // vertices from map.
1258  for (perIt = periodicVerts.begin();
1259  perIt != periodicVerts.end(); ++perIt)
1260  {
1261  if (locVerts.count(perIt->first) > 0)
1262  {
1263  m_periodicVerts.insert(*perIt);
1264  }
1265  }
1266  }
1267 
1268  bool DisContField2D::IsLeftAdjacentEdge(const int n, const int e)
1269  {
1270  set<int>::iterator it;
1272  m_traceMap->GetElmtToTrace()[n][e]->
1273  as<LocalRegions::Expansion1D>();
1274 
1275 
1276  bool fwd = true;
1277  if (traceEl->GetLeftAdjacentElementEdge () == -1 ||
1278  traceEl->GetRightAdjacentElementEdge() == -1)
1279  {
1280  // Boundary edge (1 connected element). Do nothing in
1281  // serial.
1282  it = m_boundaryEdges.find(traceEl->GetElmtId());
1283 
1284  // If the edge does not have a boundary condition set on
1285  // it, then assume it is a partition edge.
1286  if (it == m_boundaryEdges.end())
1287  {
1288  int traceGeomId = traceEl->GetGeom1D()->GetGlobalID();
1290  traceGeomId);
1291 
1292  if (pIt != m_periodicEdges.end() && !pIt->second[0].isLocal)
1293  {
1294  fwd = traceGeomId == min(traceGeomId,pIt->second[0].id);
1295  }
1296  else
1297  {
1298  int offset = m_trace->GetPhys_Offset(traceEl->GetElmtId());
1299 
1300  fwd = m_traceMap->
1301  GetTraceToUniversalMapUnique(offset) >= 0;
1302  }
1303  }
1304  }
1305  else if (traceEl->GetLeftAdjacentElementEdge () != -1 &&
1306  traceEl->GetRightAdjacentElementEdge() != -1)
1307  {
1308  // Non-boundary edge (2 connected elements).
1309  fwd = dynamic_cast<Nektar::StdRegions::StdExpansion*>
1310  (traceEl->GetLeftAdjacentElementExp().get()) ==
1311  (*m_exp)[n].get();
1312  }
1313  else
1314  {
1315  ASSERTL2(false, "Unconnected trace element!");
1316  }
1317 
1318  return fwd;
1319  }
1320 
1321  // Construct the two trace vectors of the inner and outer
1322  // trace solution from the field contained in m_phys, where
1323  // the Weak dirichlet boundary conditions are listed in the
1324  // outer part of the vecotr
1328  {
1329  v_GetFwdBwdTracePhys(m_phys, Fwd, Bwd);
1330  }
1331 
1332  /**
1333  * @brief This method extracts the "forward" and "backward" trace data
1334  * from the array @a field and puts the data into output vectors @a Fwd
1335  * and @a Bwd.
1336  *
1337  * We first define the convention which defines "forwards" and
1338  * "backwards". First an association is made between the edge of each
1339  * element and its corresponding edge in the trace space using the
1340  * mapping #m_traceMap. The element can either be left-adjacent or
1341  * right-adjacent to this trace edge (see
1342  * Expansion1D::GetLeftAdjacentElementExp). Boundary edges are always
1343  * left-adjacent since left-adjacency is populated first.
1344  *
1345  * If the element is left-adjacent we extract the edge trace data from
1346  * @a field into the forward trace space @a Fwd; otherwise, we place it
1347  * in the backwards trace space @a Bwd. In this way, we form a unique
1348  * set of trace normals since these are always extracted from
1349  * left-adjacent elements.
1350  *
1351  * @param field is a NekDouble array which contains the 2D data
1352  * from which we wish to extract the backward and
1353  * forward orientated trace/edge arrays.
1354  * @param Fwd The resulting forwards space.
1355  * @param Bwd The resulting backwards space.
1356  */
1358  const Array<OneD, const NekDouble> &field,
1361  {
1362  int cnt, n, e, npts, phys_offset;
1363 
1364  // Zero forward/backward vectors.
1365  Vmath::Zero(Fwd.num_elements(), Fwd, 1);
1366  Vmath::Zero(Bwd.num_elements(), Bwd, 1);
1367 
1368  // Basis definition on each element
1369  LibUtilities::BasisSharedPtr basis = (*m_exp)[0]->GetBasis(0);
1370  if (basis->GetBasisType() != LibUtilities::eGauss_Lagrange)
1371  {
1372 
1373  // blocked routine
1375  GetNLocTracePts());
1376 
1377  m_locTraceToTraceMap->LocTracesFromField(field, edgevals);
1378  m_locTraceToTraceMap->InterpLocEdgesToTrace(0, edgevals, Fwd);
1379 
1380  Array<OneD, NekDouble> invals = edgevals + m_locTraceToTraceMap->
1381  GetNFwdLocTracePts();
1382  m_locTraceToTraceMap->InterpLocEdgesToTrace(1, invals, Bwd);
1383  }
1384  else
1385  {
1386  // Loop over elements and collect forward expansion
1387  int nexp = GetExpSize();
1388  Array<OneD,NekDouble> e_tmp;
1390  boost::unordered_map<int,pair<int,int> >::iterator it3;
1392 
1394  &elmtToTrace = m_traceMap->GetElmtToTrace();
1395 
1396  for(cnt = n = 0; n < nexp; ++n)
1397  {
1398  exp2d = (*m_exp)[n]->as<LocalRegions::Expansion2D>();
1399  phys_offset = GetPhys_Offset(n);
1400 
1401  for(e = 0; e < exp2d->GetNedges(); ++e, ++cnt)
1402  {
1403  int offset = m_trace->GetPhys_Offset(
1404  elmtToTrace[n][e]->GetElmtId());
1405 
1406  if (m_leftAdjacentEdges[cnt])
1407  {
1408  exp2d->GetEdgePhysVals(e, elmtToTrace[n][e],
1409  field + phys_offset,
1410  e_tmp = Fwd + offset);
1411  }
1412  else
1413  {
1414  exp2d->GetEdgePhysVals(e, elmtToTrace[n][e],
1415  field + phys_offset,
1416  e_tmp = Bwd + offset);
1417  }
1418 
1419  }
1420  }
1421  }
1422 
1423  // Fill boundary conditions into missing elements
1424  int id1, id2 = 0;
1425 
1426  for (cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
1427  {
1428  if (m_bndConditions[n]->GetBoundaryConditionType() ==
1430  {
1431  for (e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
1432  {
1433  npts = m_bndCondExpansions[n]->GetExp(e)->GetNumPoints(0);
1434  id1 = m_bndCondExpansions[n]->GetPhys_Offset(e);
1435  id2 = m_trace->GetPhys_Offset(m_traceMap->
1436  GetBndCondTraceToGlobalTraceMap(cnt+e));
1437  Vmath::Vcopy(npts,
1438  &(m_bndCondExpansions[n]->GetPhys())[id1], 1,
1439  &Bwd[id2], 1);
1440  }
1441 
1442  cnt += e;
1443  }
1444  else if (m_bndConditions[n]->GetBoundaryConditionType() ==
1446  m_bndConditions[n]->GetBoundaryConditionType() ==
1448  {
1449  for(e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
1450  {
1451  npts = m_bndCondExpansions[n]->GetExp(e)->GetNumPoints(0);
1452  id1 = m_bndCondExpansions[n]->GetPhys_Offset(e);
1453  ASSERTL0((m_bndCondExpansions[n]->GetPhys())[id1]==0.0,
1454  "Method not set up for non-zero Neumann "
1455  "boundary condition");
1456  id2 = m_trace->GetPhys_Offset(
1457  m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
1458  Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
1459  }
1460 
1461  cnt += e;
1462  }
1463  else if (m_bndConditions[n]->GetBoundaryConditionType() !=
1465  {
1466  ASSERTL0(false,
1467  "Method not set up for this boundary condition.");
1468  }
1469  }
1470 
1471  // Copy any periodic boundary conditions.
1472  for (n = 0; n < m_periodicFwdCopy.size(); ++n)
1473  {
1474  Bwd[m_periodicBwdCopy[n]] = Fwd[m_periodicFwdCopy[n]];
1475  }
1476 
1477  // Do parallel exchange for forwards/backwards spaces.
1478  m_traceMap->UniversalTraceAssemble(Fwd);
1479  m_traceMap->UniversalTraceAssemble(Bwd);
1480  }
1481 
1482 
1484  Array<OneD, NekDouble> &outarray)
1485  {
1486  ASSERTL1(m_physState == true,
1487  "Field must be in physical state to extract trace space.");
1488 
1489  v_ExtractTracePhys(m_phys, outarray);
1490  }
1491 
1492  /**
1493  * @brief This method extracts the trace (edges in 2D) from the field @a
1494  * inarray and puts the values in @a outarray.
1495  *
1496  * It assumes the field is C0 continuous so that it can overwrite the
1497  * edge data when visited by the two adjacent elements.
1498  *
1499  * @param inarray An array containing the 2D data from which we wish
1500  * to extract the edge data.
1501  * @param outarray The resulting edge information.
1502  */
1504  const Array<OneD, const NekDouble> &inarray,
1505  Array<OneD, NekDouble> &outarray)
1506  {
1507  // Loop over elemente and collect forward expansion
1508  int nexp = GetExpSize();
1509  int n, e, offset, phys_offset;
1510  Array<OneD,NekDouble> e_tmp;
1512  &elmtToTrace = m_traceMap->GetElmtToTrace();
1513 
1514  ASSERTL1(outarray.num_elements() >= m_trace->GetNpoints(),
1515  "input array is of insufficient length");
1516 
1517  // use m_trace tmp space in element to fill values
1518  for(n = 0; n < nexp; ++n)
1519  {
1520  phys_offset = GetPhys_Offset(n);
1521 
1522  for(e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
1523  {
1524  offset = m_trace->GetPhys_Offset(
1525  elmtToTrace[n][e]->GetElmtId());
1526  (*m_exp)[n]->GetEdgePhysVals(e, elmtToTrace[n][e],
1527  inarray + phys_offset,
1528  e_tmp = outarray + offset);
1529  }
1530  }
1531  }
1532 
1534  const Array<OneD, const NekDouble> &Fx,
1535  const Array<OneD, const NekDouble> &Fy,
1536  Array<OneD, NekDouble> &outarray)
1537  {
1538  int e, n, offset, t_offset;
1539  Array<OneD, NekDouble> e_outarray;
1541  &elmtToTrace = m_traceMap->GetElmtToTrace();
1542 
1543  for(n = 0; n < GetExpSize(); ++n)
1544  {
1545  offset = GetCoeff_Offset(n);
1546  for(e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
1547  {
1548  t_offset = GetTrace()->GetPhys_Offset(
1549  elmtToTrace[n][e]->GetElmtId());
1550 
1551  (*m_exp)[n]->AddEdgeNormBoundaryInt(
1552  e,elmtToTrace[n][e],
1553  Fx + t_offset,
1554  Fy + t_offset,
1555  e_outarray = outarray+offset);
1556  }
1557  }
1558  }
1559 
1560  /**
1561  * @brief Add trace contributions into elemental coefficient spaces.
1562  *
1563  * Given some quantity \f$ \vec{Fn} \f$, which conatins this
1564  * routine calculates the integral
1565  *
1566  * \f[
1567  * \int_{\Omega^e} \vec{Fn}, \mathrm{d}S
1568  * \f]
1569  *
1570  * and adds this to the coefficient space provided by outarray.
1571  *
1572  * @see Expansion2D::AddEdgeNormBoundaryInt
1573  *
1574  * @param Fn The trace quantities.
1575  * @param outarray Resulting 2D coefficient space.
1576  */
1578  const Array<OneD, const NekDouble> &Fn,
1579  Array<OneD, NekDouble> &outarray)
1580  {
1581  // Basis definition on each element
1582  LibUtilities::BasisSharedPtr basis = (*m_exp)[0]->GetBasis(0);
1583  if (basis->GetBasisType() != LibUtilities::eGauss_Lagrange)
1584  {
1585  Array<OneD, NekDouble> Fcoeffs(m_trace->GetNcoeffs());
1586  m_trace->IProductWRTBase(Fn, Fcoeffs);
1587 
1588  m_locTraceToTraceMap->AddTraceCoeffsToFieldCoeffs(Fcoeffs,
1589  outarray);
1590  }
1591  else
1592  {
1593  int e, n, offset, t_offset;
1594  Array<OneD, NekDouble> e_outarray;
1596  &elmtToTrace = m_traceMap->GetElmtToTrace();
1597 
1598  for(n = 0; n < GetExpSize(); ++n)
1599  {
1600  offset = GetCoeff_Offset(n);
1601  for(e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
1602  {
1603  t_offset = GetTrace()->GetPhys_Offset(
1604  elmtToTrace[n][e]->GetElmtId());
1605  (*m_exp)[n]->AddEdgeNormBoundaryInt(
1606  e, elmtToTrace[n][e],
1607  Fn+t_offset,
1608  e_outarray = outarray+offset);
1609  }
1610  }
1611  }
1612  }
1613 
1614 
1615  /**
1616  * @brief Add trace contributions into elemental coefficient spaces.
1617  *
1618  * Given some quantity \f$ \vec{q} \f$, calculate the elemental integral
1619  *
1620  * \f[
1621  * \int_{\Omega^e} \vec{q}, \mathrm{d}S
1622  * \f]
1623  *
1624  * and adds this to the coefficient space provided by
1625  * outarray. The value of q is determined from the routine
1626  * IsLeftAdjacentEdge() which if true we use Fwd else we use
1627  * Bwd
1628  *
1629  * @see Expansion2D::AddEdgeNormBoundaryInt
1630  *
1631  * @param Fwd The trace quantities associated with left (fwd)
1632  * adjancent elmt.
1633  * @param Bwd The trace quantities associated with right (bwd)
1634  * adjacent elet.
1635  * @param outarray Resulting 2D coefficient space.
1636  */
1638  const Array<OneD, const NekDouble> &Fwd,
1639  const Array<OneD, const NekDouble> &Bwd,
1640  Array<OneD, NekDouble> &outarray)
1641  {
1642  int e,n,offset, t_offset;
1643  Array<OneD, NekDouble> e_outarray;
1645  &elmtToTrace = m_traceMap->GetElmtToTrace();
1646 
1647  for (n = 0; n < GetExpSize(); ++n)
1648  {
1649  offset = GetCoeff_Offset(n);
1650  for (e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
1651  {
1652  t_offset = GetTrace()->GetPhys_Offset(
1653  elmtToTrace[n][e]->GetElmtId());
1654 
1655  // Evaluate upwind flux less local edge
1656  if (IsLeftAdjacentEdge(n, e))
1657  {
1658  (*m_exp)[n]->AddEdgeNormBoundaryInt(
1659  e, elmtToTrace[n][e], Fwd+t_offset,
1660  e_outarray = outarray+offset);
1661  }
1662  else
1663  {
1664  (*m_exp)[n]->AddEdgeNormBoundaryInt(
1665  e, elmtToTrace[n][e], Bwd+t_offset,
1666  e_outarray = outarray+offset);
1667  }
1668 
1669  }
1670  }
1671  }
1672 
1673  /**
1674  * @brief Set up a list of element IDs and edge IDs that link to the
1675  * boundary conditions.
1676  */
1678  Array<OneD, int> &ElmtID,
1679  Array<OneD, int> &EdgeID)
1680  {
1681  if (m_BCtoElmMap.num_elements() == 0)
1682  {
1683  map<int, int> globalIdMap;
1684  int i,n;
1685  int cnt;
1686  int nbcs = 0;
1687 
1689  boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(
1690  m_graph);
1691 
1692  // Populate global ID map (takes global geometry ID to local
1693  // expansion list ID).
1694  for (i = 0; i < GetExpSize(); ++i)
1695  {
1696  globalIdMap[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
1697  }
1698 
1699  // Determine number of boundary condition expansions.
1700  for(i = 0; i < m_bndConditions.num_elements(); ++i)
1701  {
1702  nbcs += m_bndCondExpansions[i]->GetExpSize();
1703  }
1704 
1705  // Initialize arrays
1708 
1710  for (cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
1711  {
1712  for (i = 0; i < m_bndCondExpansions[n]->GetExpSize();
1713  ++i, ++cnt)
1714  {
1715  exp1d = m_bndCondExpansions[n]->GetExp(i)->
1716  as<LocalRegions::Expansion1D>();
1717  // Use edge to element map from MeshGraph2D.
1719  graph2D->GetElementsFromEdge(exp1d->GetGeom1D());
1720 
1721  m_BCtoElmMap[cnt] = globalIdMap[(*tmp)[0]->
1722  m_Element->GetGlobalID()];
1723  m_BCtoEdgMap[cnt] = (*tmp)[0]->m_EdgeIndx;
1724  }
1725  }
1726  }
1727  ElmtID = m_BCtoElmMap;
1728  EdgeID = m_BCtoEdgMap;
1729  }
1730 
1732  boost::shared_ptr<ExpList> &result,
1733  const bool DeclareCoeffPhysArrays)
1734  {
1735  int n, cnt, nq;
1736  int offsetOld, offsetNew;
1737  std::vector<unsigned int> eIDs;
1738 
1739  Array<OneD, int> ElmtID,EdgeID;
1740  GetBoundaryToElmtMap(ElmtID,EdgeID);
1741 
1742  // Skip other boundary regions
1743  for (cnt = n = 0; n < i; ++n)
1744  {
1745  cnt += m_bndCondExpansions[n]->GetExpSize();
1746  }
1747 
1748  // Populate eIDs with information from BoundaryToElmtMap
1749  for (n = 0; n < m_bndCondExpansions[i]->GetExpSize(); ++n)
1750  {
1751  eIDs.push_back(ElmtID[cnt+n]);
1752  }
1753 
1754  // Create expansion list
1755  result =
1757  (*this, eIDs, DeclareCoeffPhysArrays);
1758 
1759  // Copy phys and coeffs to new explist
1760  if( DeclareCoeffPhysArrays)
1761  {
1762  Array<OneD, NekDouble> tmp1, tmp2;
1763  for (n = 0; n < result->GetExpSize(); ++n)
1764  {
1765  nq = GetExp(ElmtID[cnt+n])->GetTotPoints();
1766  offsetOld = GetPhys_Offset(ElmtID[cnt+n]);
1767  offsetNew = result->GetPhys_Offset(n);
1768  Vmath::Vcopy(nq, tmp1 = GetPhys()+ offsetOld, 1,
1769  tmp2 = result->UpdatePhys()+ offsetNew, 1);
1770 
1771  nq = GetExp(ElmtID[cnt+n])->GetNcoeffs();
1772  offsetOld = GetCoeff_Offset(ElmtID[cnt+n]);
1773  offsetNew = result->GetCoeff_Offset(n);
1774  Vmath::Vcopy(nq, tmp1 = GetCoeffs()+ offsetOld, 1,
1775  tmp2 = result->UpdateCoeffs()+ offsetNew, 1);
1776  }
1777  }
1778  }
1779 
1780  /**
1781  * @brief Reset this field, so that geometry information can be updated.
1782  */
1784  {
1785  ExpList::v_Reset();
1786 
1787  // Reset boundary condition expansions.
1788  for (int n = 0; n < m_bndCondExpansions.num_elements(); ++n)
1789  {
1790  m_bndCondExpansions[n]->Reset();
1791  }
1792  }
1793 
1794  /**
1795  * @brief Calculate the \f$ L^2 \f$ error of the \f$ Q_{\rm dir} \f$
1796  * derivative using the consistent DG evaluation of \f$ Q_{\rm dir} \f$.
1797  *
1798  * The solution provided is of the primative variation at the quadrature
1799  * points and the derivative is compared to the discrete derivative at
1800  * these points, which is likely to be undesirable unless using a much
1801  * higher number of quadrature points than the polynomial order used to
1802  * evaluate \f$ Q_{\rm dir} \f$.
1803  */
1805  const int dir,
1806  const Array<OneD, const NekDouble> &soln)
1807  {
1808  int i,e,ncoeff_edge;
1809  Array<OneD, const NekDouble> tmp_coeffs;
1810  Array<OneD, NekDouble> out_d(m_ncoeffs), out_tmp;
1811 
1813  &elmtToTrace = m_traceMap->GetElmtToTrace();
1814 
1815  StdRegions::Orientation edgedir;
1816 
1817  int eid,cnt;
1818  int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
1819  Array<OneD, NekDouble> loc_lambda(LocBndCoeffs), edge_lambda;
1820  m_traceMap->GlobalToLocalBnd(m_trace->GetCoeffs(),loc_lambda);
1821 
1822  edge_lambda = loc_lambda;
1823 
1824  // Calculate Q using standard DG formulation.
1825  for(i = cnt = 0; i < GetExpSize(); ++i)
1826  {
1827  eid = m_offset_elmt_id[i];
1828 
1829  // Probably a better way of setting up lambda than this.
1830  // Note cannot use PutCoeffsInToElmts since lambda space
1831  // is mapped during the solve.
1832  int nEdges = (*m_exp)[i]->GetNedges();
1833  Array<OneD, Array<OneD, NekDouble> > edgeCoeffs(nEdges);
1834 
1835  for(e = 0; e < nEdges; ++e)
1836  {
1837  edgedir = (*m_exp)[eid]->GetEorient(e);
1838  ncoeff_edge = elmtToTrace[eid][e]->GetNcoeffs();
1839  edgeCoeffs[e] = Array<OneD, NekDouble>(ncoeff_edge);
1840  Vmath::Vcopy(ncoeff_edge, edge_lambda, 1, edgeCoeffs[e], 1);
1841  elmtToTrace[eid][e]->SetCoeffsToOrientation(
1842  edgedir, edgeCoeffs[e], edgeCoeffs[e]);
1843  edge_lambda = edge_lambda + ncoeff_edge;
1844  }
1845 
1846  (*m_exp)[eid]->DGDeriv(dir,
1847  tmp_coeffs=m_coeffs+m_coeff_offset[eid],
1848  elmtToTrace[eid],
1849  edgeCoeffs,
1850  out_tmp = out_d+cnt);
1851  cnt += (*m_exp)[eid]->GetNcoeffs();
1852  }
1853 
1854  BwdTrans(out_d,m_phys);
1855  Vmath::Vsub(m_npoints,m_phys,1,soln,1,m_phys,1);
1856  return L2(m_phys);
1857  }
1858 
1860  const Array<OneD, const NekDouble> &inarray,
1861  Array<OneD, NekDouble> &outarray,
1862  const FlagList &flags,
1863  const StdRegions::ConstFactorMap &factors,
1864  const StdRegions::VarCoeffMap &varcoeff,
1865  const Array<OneD, const NekDouble> &dirForcing,
1866  const bool PhysSpaceForcing)
1867 
1868  {
1869  int i,j,n,cnt,cnt1,nbndry;
1870  int nexp = GetExpSize();
1872 
1874  DNekVec F(m_ncoeffs,f,eWrapper);
1875  Array<OneD,NekDouble> e_f, e_l;
1876 
1877  //----------------------------------
1878  // Setup RHS Inner product
1879  //----------------------------------
1880  if(PhysSpaceForcing)
1881  {
1882  IProductWRTBase(inarray,f);
1883  Vmath::Neg(m_ncoeffs,f,1);
1884  }
1885  else
1886  {
1887  Vmath::Smul(m_ncoeffs,-1.0,inarray,1,f,1);
1888  }
1889 
1890  //----------------------------------
1891  // Solve continuous flux System
1892  //----------------------------------
1893  int GloBndDofs = m_traceMap->GetNumGlobalBndCoeffs();
1894  int NumDirichlet = m_traceMap->GetNumLocalDirBndCoeffs();
1895  int e_ncoeffs,id;
1896 
1897  // Retrieve block matrix of U^e
1899  NullAssemblyMapSharedPtr,factors,varcoeff);
1900  const DNekScalBlkMatSharedPtr &HDGLamToU = GetBlockMatrix(
1901  HDGLamToUKey);
1902 
1903  // Retrieve global trace space storage, \Lambda, from trace
1904  // expansion
1905  Array<OneD,NekDouble> BndSol = m_trace->UpdateCoeffs();
1906 
1907  // Create trace space forcing, F
1908  Array<OneD,NekDouble> BndRhs(GloBndDofs,0.0);
1909 
1910  // Zero \Lambda
1911  Vmath::Zero(GloBndDofs,BndSol,1);
1912 
1913  // Retrieve number of local trace space coefficients N_{\lambda},
1914  // and set up local elemental trace solution \lambda^e.
1915  int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
1916  Array<OneD, NekDouble> loc_lambda(LocBndCoeffs);
1917  DNekVec LocLambda(LocBndCoeffs,loc_lambda,eWrapper);
1918 
1919  //----------------------------------
1920  // Evaluate Trace Forcing vector F
1921  // Kirby et al, 2010, P23, Step 5.
1922  //----------------------------------
1923  // Loop over all expansions in the domain
1924  for(cnt = cnt1 = n = 0; n < nexp; ++n)
1925  {
1926  nbndry = (*m_exp)[m_offset_elmt_id[n]]->NumDGBndryCoeffs();
1927 
1928  e_ncoeffs = (*m_exp)[m_offset_elmt_id[n]]->GetNcoeffs();
1929  e_f = f + cnt;
1930  e_l = loc_lambda + cnt1;
1931 
1932  // Local trace space \lambda^e
1933  DNekVec Floc (nbndry, e_l, eWrapper);
1934  // Local forcing f^e
1935  DNekVec ElmtFce (e_ncoeffs, e_f, eWrapper);
1936  // Compute local (U^e)^{\top} f^e
1937  Floc = Transpose(*(HDGLamToU->GetBlock(n,n)))*ElmtFce;
1938 
1939  cnt += e_ncoeffs;
1940  cnt1 += nbndry;
1941  }
1942 
1943  // Assemble local \lambda_e into global \Lambda
1944  m_traceMap->AssembleBnd(loc_lambda,BndRhs);
1945 
1946  // Copy Dirichlet boundary conditions and weak forcing into trace
1947  // space
1948  cnt = 0;
1949  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
1950  {
1951  if(m_bndConditions[i]->GetBoundaryConditionType() ==
1953  {
1954  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
1955  {
1956  id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(cnt++);
1957  BndSol[id] = m_bndCondExpansions[i]->GetCoeffs()[j];
1958  }
1959  }
1960  else
1961  {
1962  //Add weak boundary condition to trace forcing
1963  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
1964  {
1965  id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(cnt++);
1966  BndRhs[id] += m_bndCondExpansions[i]->GetCoeffs()[j];
1967  }
1968  }
1969  }
1970 
1971  //----------------------------------
1972  // Solve trace problem: \Lambda = K^{-1} F
1973  // K is the HybridDGHelmBndLam matrix.
1974  //----------------------------------
1975  if(GloBndDofs - NumDirichlet > 0)
1976  {
1978  m_traceMap,factors,varcoeff);
1980  LinSys->Solve(BndRhs,BndSol,m_traceMap);
1981  }
1982 
1983  //----------------------------------
1984  // Internal element solves
1985  //----------------------------------
1987  NullAssemblyMapSharedPtr,factors,varcoeff);
1988  const DNekScalBlkMatSharedPtr& InvHDGHelm = GetBlockMatrix(
1989  invHDGhelmkey);
1990  DNekVec out(m_ncoeffs,outarray,eWrapper);
1991  Vmath::Zero(m_ncoeffs,outarray,1);
1992 
1993  // get local trace solution from BndSol
1994  m_traceMap->GlobalToLocalBnd(BndSol,loc_lambda);
1995 
1996  // out = u_f + u_lam = (*InvHDGHelm)*f + (LamtoU)*Lam
1997  out = (*InvHDGHelm)*F + (*HDGLamToU)*LocLambda;
1998  }
1999 
2000 
2001  /**
2002  * @brief Calculates the result of the multiplication of a global matrix
2003  * of type specified by @a mkey with a vector given by @a inarray.
2004  *
2005  * @param mkey Key representing desired matrix multiplication.
2006  * @param inarray Input vector.
2007  * @param outarray Resulting multiplication.
2008  */
2010  const GlobalMatrixKey &gkey,
2011  const Array<OneD,const NekDouble> &inarray,
2012  Array<OneD, NekDouble> &outarray,
2013  CoeffState coeffstate)
2014  {
2015  int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
2016  Array<OneD, NekDouble> loc_lambda(LocBndCoeffs);
2017  DNekVec LocLambda(LocBndCoeffs,loc_lambda,eWrapper);
2018  const DNekScalBlkMatSharedPtr& HDGHelm = GetBlockMatrix(gkey);
2019 
2020  m_traceMap->GlobalToLocalBnd(inarray, loc_lambda);
2021  LocLambda = (*HDGHelm) * LocLambda;
2022  m_traceMap->AssembleBnd(loc_lambda,outarray);
2023  }
2024 
2025  /**
2026  * @brief Search through the edge expansions and identify which ones
2027  * have Robin/Mixed type boundary conditions.
2028  *
2029  * If a Robin boundary is found then store the edge ID of the boundary
2030  * condition and the array of points of the physical space boundary
2031  * condition which are hold the boundary condition primitive variable
2032  * coefficient at the quatrature points
2033  *
2034  * @return A map containing the Robin boundary condition information
2035  * using a key of the element ID.
2036  */
2037  map<int, RobinBCInfoSharedPtr> DisContField2D::v_GetRobinBCInfo(void)
2038  {
2039  int i,cnt;
2040  map<int, RobinBCInfoSharedPtr> returnval;
2041  Array<OneD, int> ElmtID,EdgeID;
2042  GetBoundaryToElmtMap(ElmtID,EdgeID);
2043 
2044  for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
2045  {
2046  MultiRegions::ExpListSharedPtr locExpList;
2047 
2048  if(m_bndConditions[i]->GetBoundaryConditionType() ==
2050  {
2051  int e,elmtid;
2052  Array<OneD, NekDouble> Array_tmp;
2053 
2054  locExpList = m_bndCondExpansions[i];
2055 
2056  int npoints = locExpList->GetNpoints();
2057  Array<OneD, NekDouble> x0(npoints, 0.0);
2058  Array<OneD, NekDouble> x1(npoints, 0.0);
2059  Array<OneD, NekDouble> x2(npoints, 0.0);
2060  Array<OneD, NekDouble> coeffphys(npoints);
2061 
2062  locExpList->GetCoords(x0, x1, x2);
2063 
2064  LibUtilities::Equation coeffeqn =
2065  boost::static_pointer_cast<
2067  (m_bndConditions[i])->m_robinPrimitiveCoeff;
2068 
2069  // evalaute coefficient
2070  coeffeqn.Evaluate(x0, x1, x2, 0.0, coeffphys);
2071 
2072  for(e = 0; e < locExpList->GetExpSize(); ++e)
2073  {
2074  RobinBCInfoSharedPtr rInfo =
2077  EdgeID[cnt+e],
2078  Array_tmp = coeffphys +
2079  locExpList->GetPhys_Offset(e));
2080 
2081  elmtid = ElmtID[cnt+e];
2082  // make link list if necessary
2083  if(returnval.count(elmtid) != 0)
2084  {
2085  rInfo->next = returnval.find(elmtid)->second;
2086  }
2087  returnval[elmtid] = rInfo;
2088  }
2089  }
2090  cnt += m_bndCondExpansions[i]->GetExpSize();
2091  }
2092 
2093  return returnval;
2094  }
2095 
2096  /**
2097  * @brief Evaluate HDG post-processing to increase polynomial order of
2098  * solution.
2099  *
2100  * This function takes the solution (assumed to be one order lower) in
2101  * physical space, and postprocesses at the current polynomial order by
2102  * solving the system:
2103  *
2104  * \f[
2105  * \begin{aligned}
2106  * (\nabla w, \nabla u^*) &= (\nabla w, u), \\
2107  * \langle \nabla u^*, 1 \rangle &= \langle \nabla u, 1 \rangle
2108  * \end{aligned}
2109  * \f]
2110  *
2111  * where \f$ u \f$ corresponds with the current solution as stored
2112  * inside #m_coeffs.
2113  *
2114  * @param outarray The resulting field \f$ u^* \f$.
2115  */
2117  Array<OneD, NekDouble> &outarray)
2118  {
2119  int i,cnt,e,ncoeff_edge;
2120  Array<OneD, NekDouble> force, out_tmp, qrhs, qrhs1;
2122  &elmtToTrace = m_traceMap->GetElmtToTrace();
2123 
2124  StdRegions::Orientation edgedir;
2125 
2126  int eid,nq_elmt, nm_elmt;
2127  int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
2128  Array<OneD, NekDouble> loc_lambda(LocBndCoeffs), edge_lambda;
2129  Array<OneD, NekDouble> tmp_coeffs;
2130  m_traceMap->GlobalToLocalBnd(m_trace->GetCoeffs(),loc_lambda);
2131 
2132  edge_lambda = loc_lambda;
2133 
2134  // Calculate Q using standard DG formulation.
2135  for(i = cnt = 0; i < GetExpSize(); ++i)
2136  {
2137  eid = m_offset_elmt_id[i];
2138 
2139  nq_elmt = (*m_exp)[eid]->GetTotPoints();
2140  nm_elmt = (*m_exp)[eid]->GetNcoeffs();
2141  qrhs = Array<OneD, NekDouble>(nq_elmt);
2142  qrhs1 = Array<OneD, NekDouble>(nq_elmt);
2143  force = Array<OneD, NekDouble>(2*nm_elmt);
2144  out_tmp = force + nm_elmt;
2146 
2147  int num_points0 = (*m_exp)[eid]->GetBasis(0)->GetNumPoints();
2148  int num_points1 = (*m_exp)[eid]->GetBasis(1)->GetNumPoints();
2149  int num_modes0 = (*m_exp)[eid]->GetBasis(0)->GetNumModes();
2150  int num_modes1 = (*m_exp)[eid]->GetBasis(1)->GetNumModes();
2151 
2152  // Probably a better way of setting up lambda than this. Note
2153  // cannot use PutCoeffsInToElmts since lambda space is mapped
2154  // during the solve.
2155  int nEdges = (*m_exp)[i]->GetNedges();
2156  Array<OneD, Array<OneD, NekDouble> > edgeCoeffs(nEdges);
2157 
2158  for(e = 0; e < (*m_exp)[eid]->GetNedges(); ++e)
2159  {
2160  edgedir = (*m_exp)[eid]->GetEorient(e);
2161  ncoeff_edge = elmtToTrace[eid][e]->GetNcoeffs();
2162  edgeCoeffs[e] = Array<OneD, NekDouble>(ncoeff_edge);
2163  Vmath::Vcopy(ncoeff_edge, edge_lambda, 1, edgeCoeffs[e], 1);
2164  elmtToTrace[eid][e]->SetCoeffsToOrientation(
2165  edgedir, edgeCoeffs[e], edgeCoeffs[e]);
2166  edge_lambda = edge_lambda + ncoeff_edge;
2167  }
2168 
2169  //creating orthogonal expansion (checking if we have quads or triangles)
2170  LibUtilities::ShapeType shape = (*m_exp)[eid]->DetShapeType();
2171  switch(shape)
2172  {
2174  {
2177  LibUtilities::BasisKey BkeyQ1(LibUtilities::eOrtho_A, num_modes0, PkeyQ1);
2178  LibUtilities::BasisKey BkeyQ2(LibUtilities::eOrtho_A, num_modes1, PkeyQ2);
2179  SpatialDomains::QuadGeomSharedPtr qGeom = boost::dynamic_pointer_cast<SpatialDomains::QuadGeom>((*m_exp)[eid]->GetGeom());
2180  ppExp = MemoryManager<LocalRegions::QuadExp>::AllocateSharedPtr(BkeyQ1, BkeyQ2, qGeom);
2181  }
2182  break;
2184  {
2187  LibUtilities::BasisKey BkeyT1(LibUtilities::eOrtho_A, num_modes0, PkeyT1);
2188  LibUtilities::BasisKey BkeyT2(LibUtilities::eOrtho_B, num_modes1, PkeyT2);
2189  SpatialDomains::TriGeomSharedPtr tGeom = boost::dynamic_pointer_cast<SpatialDomains::TriGeom>((*m_exp)[eid]->GetGeom());
2190  ppExp = MemoryManager<LocalRegions::TriExp>::AllocateSharedPtr(BkeyT1, BkeyT2, tGeom);
2191  }
2192  break;
2193  default:
2194  ASSERTL0(false, "Wrong shape type, HDG postprocessing is not implemented");
2195  };
2196 
2197 
2198  //DGDeriv
2199  // (d/dx w, d/dx q_0)
2200  (*m_exp)[eid]->DGDeriv(
2201  0,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
2202  elmtToTrace[eid], edgeCoeffs, out_tmp);
2203  (*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
2204  //(*m_exp)[eid]->IProductWRTDerivBase(0,qrhs,force);
2205  ppExp->IProductWRTDerivBase(0,qrhs,force);
2206 
2207 
2208  // + (d/dy w, d/dy q_1)
2209  (*m_exp)[eid]->DGDeriv(
2210  1,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
2211  elmtToTrace[eid], edgeCoeffs, out_tmp);
2212 
2213  (*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
2214  //(*m_exp)[eid]->IProductWRTDerivBase(1,qrhs,out_tmp);
2215  ppExp->IProductWRTDerivBase(1,qrhs,out_tmp);
2216 
2217  Vmath::Vadd(nm_elmt,force,1,out_tmp,1,force,1);
2218 
2219  // determine force[0] = (1,u)
2220  (*m_exp)[eid]->BwdTrans(
2221  tmp_coeffs = m_coeffs + m_coeff_offset[eid],qrhs);
2222  force[0] = (*m_exp)[eid]->Integral(qrhs);
2223 
2224  // multiply by inverse Laplacian matrix
2225  // get matrix inverse
2226  LocalRegions::MatrixKey lapkey(StdRegions::eInvLaplacianWithUnityMean, ppExp->DetShapeType(), *ppExp);
2227  DNekScalMatSharedPtr lapsys = ppExp->GetLocMatrix(lapkey);
2228 
2229  NekVector<NekDouble> in (nm_elmt,force,eWrapper);
2230  NekVector<NekDouble> out(nm_elmt);
2231 
2232  out = (*lapsys)*in;
2233 
2234  // Transforming back to modified basis
2235  Array<OneD, NekDouble> work(nq_elmt);
2236  ppExp->BwdTrans(out.GetPtr(), work);
2237  (*m_exp)[eid]->FwdTrans(work, tmp_coeffs = outarray + m_coeff_offset[eid]);
2238  }
2239  }
2240 
2241  /**
2242  * Evaluates the boundary condition expansions, \a bndCondExpansions,
2243  * given the information provided by \a bndConditions.
2244  * @param time The time at which the boundary conditions
2245  * should be evaluated.
2246  * @param bndCondExpansions List of boundary conditions.
2247  * @param bndConditions Information about the boundary conditions.
2248  *
2249  * This will only be undertaken for time dependent
2250  * boundary conditions unless time == 0.0 which is the
2251  * case when the method is called from the constructor.
2252  */
2254  const NekDouble time,
2255  const std::string varName,
2256  const NekDouble x2_in,
2257  const NekDouble x3_in)
2258  {
2259  int i;
2260  int npoints;
2261  int nbnd = m_bndCondExpansions.num_elements();
2262 
2263  MultiRegions::ExpListSharedPtr locExpList;
2264 
2265  for (i = 0; i < nbnd; ++i)
2266  {
2267  if (time == 0.0 ||
2268  m_bndConditions[i]->IsTimeDependent())
2269  {
2270  locExpList = m_bndCondExpansions[i];
2271  npoints = locExpList->GetNpoints();
2272  Array<OneD, NekDouble> x0(npoints, 0.0);
2273  Array<OneD, NekDouble> x1(npoints, 0.0);
2274  Array<OneD, NekDouble> x2(npoints, 0.0);
2275 
2276  // Homogeneous input case for x2.
2277  if (x2_in == NekConstants::kNekUnsetDouble)
2278  {
2279  locExpList->GetCoords(x0, x1, x2);
2280  }
2281  else
2282  {
2283  locExpList->GetCoords(x0, x1, x2);
2284  Vmath::Fill(npoints, x2_in, x2, 1);
2285  }
2286 
2287  if (m_bndConditions[i]->GetBoundaryConditionType()
2289  {
2291  = boost::static_pointer_cast<
2293  m_bndConditions[i]);
2294  string filebcs = bcPtr->m_filename;
2295 
2296  if (filebcs != "")
2297  {
2298  ExtractFileBCs(filebcs, bcPtr->m_comm, varName, locExpList);
2299  }
2300  else
2301  {
2302  LibUtilities::Equation condition =
2303  boost::static_pointer_cast<
2305  (m_bndConditions[i])->
2306  m_dirichletCondition;
2307 
2308  condition.Evaluate(x0, x1, x2, time,
2309  locExpList->UpdatePhys());
2310  }
2311 
2312  locExpList->FwdTrans_BndConstrained(
2313  locExpList->GetPhys(),
2314  locExpList->UpdateCoeffs());
2315  }
2316  else if (m_bndConditions[i]->GetBoundaryConditionType()
2318  {
2319  SpatialDomains::NeumannBCShPtr bcPtr = boost::static_pointer_cast<
2321  m_bndConditions[i]);
2322  string filebcs = bcPtr->m_filename;
2323  if (filebcs != "")
2324  {
2325  ExtractFileBCs(filebcs, bcPtr->m_comm, varName, locExpList);
2326  }
2327  else
2328  {
2329  LibUtilities::Equation condition =
2330  boost::static_pointer_cast<
2332  (m_bndConditions[i])->
2333  m_neumannCondition;
2334  condition.Evaluate(x0, x1, x2, time,
2335  locExpList->UpdatePhys());
2336  }
2337 
2338  locExpList->IProductWRTBase(
2339  locExpList->GetPhys(),
2340  locExpList->UpdateCoeffs());
2341  }
2342  else if (m_bndConditions[i]->GetBoundaryConditionType()
2344  {
2345  SpatialDomains::RobinBCShPtr bcPtr = boost::static_pointer_cast<
2347  (m_bndConditions[i]);
2348  string filebcs = bcPtr->m_filename;
2349 
2350  if (filebcs != "")
2351  {
2352  ExtractFileBCs(filebcs, bcPtr->m_comm, varName, locExpList);
2353  }
2354  else
2355  {
2356  LibUtilities::Equation condition =
2357  boost::static_pointer_cast<
2359  (m_bndConditions[i])->
2360  m_robinFunction;
2361  condition.Evaluate(x0, x1, x2, time,
2362  locExpList->UpdatePhys());
2363  }
2364 
2365  locExpList->IProductWRTBase(
2366  locExpList->GetPhys(),
2367  locExpList->UpdateCoeffs());
2368  }
2369  else
2370  {
2371  ASSERTL0(false, "This type of BC not implemented yet");
2372  }
2373  }
2374  else if (boost::iequals(m_bndConditions[i]->GetUserDefined(),
2375  "MovingBody"))
2376  {
2377  locExpList = m_bndCondExpansions[i];
2378  if (m_bndConditions[i]->GetBoundaryConditionType()
2380  {
2381  locExpList->FwdTrans_IterPerExp(
2382  locExpList->GetPhys(),
2383  locExpList->UpdateCoeffs());
2384  }
2385  else
2386  {
2387  ASSERTL0(false, "This type of BC not implemented yet");
2388  }
2389  }
2390  }
2391  }
2392  } // end of namespace
2393 } //end of namespace
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
virtual void v_ExtractTracePhys(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This method extracts the trace (edges in 2D) from the field inarray and puts the values in outarray...
#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
void FindPeriodicEdges(const SpatialDomains::BoundaryConditions &bcs, const std::string &variable)
Determine the periodic edges and vertices for the given graph.
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
Definition: ExpList.cpp:1987
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
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
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.
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
std::vector< int > m_periodicFwdCopy
A vector indicating degress of freedom which need to be copied from forwards to backwards space in ca...
const BoundaryConditionCollection & GetBoundaryConditions(void) const
Definition: Conditions.h:248
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
boost::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2075
GlobalLinSysSharedPtr GetGlobalBndLinSys(const GlobalLinSysKey &mkey)
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
STL namespace.
PeriodicMap m_periodicEdges
A map which identifies pairs of periodic edges.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1699
std::map< int, std::vector< unsigned int > > BndRegionOrdering
Definition: MeshPartition.h:54
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
PeriodicMap m_periodicVerts
A map which identifies groups of periodic vertices.
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
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Add trace contributions into elemental coefficient spaces.
std::map< int, std::vector< unsigned int > > CompositeOrdering
Definition: MeshPartition.h:53
boost::shared_ptr< NeumannBoundaryCondition > NeumannBCShPtr
Definition: Conditions.h:221
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
boost::shared_ptr< MeshGraph2D > MeshGraph2DSharedPtr
Definition: MeshGraph2D.h:238
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
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo()
Search through the edge expansions and identify which ones have Robin/Mixed type boundary conditions...
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
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
virtual void v_GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Principle Orthogonal Functions .
Definition: BasisType.h:47
boost::shared_ptr< DirichletBoundaryCondition > DirichletBCShPtr
Definition: Conditions.h:220
bool SameTypeOfBoundaryConditions(const DisContField2D &In)
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1047
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
DisContField2D()
Default constructor.
The base class for all shapes.
Definition: StdExpansion.h:69
boost::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2159
LocTraceToTraceMapSharedPtr m_locTraceToTraceMap
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
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
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)
boost::shared_ptr< ExpList1D > ExpList1DSharedPtr
Shared pointer to an ExpList1D object.
Definition: ExpList1D.h:50
virtual ~DisContField2D()
Default destructor.
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:227
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
std::set< int > m_boundaryEdges
A set storing the global IDs of any boundary edges.
int id
Geometry ID of entity.
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1024
NekDouble Evaluate() const
Definition: Equation.h:102
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
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)
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
Definition: ExpList.h:1058
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.
static std::string npts
Definition: InputFld.cpp:43
Principle Orthogonal Functions .
Definition: BasisType.h:46
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:972
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:969
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
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
bool IsLeftAdjacentEdge(const int n, const int e)
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
virtual void v_Reset()
Reset this field, so that geometry information can be updated.
static const NekDouble kNekUnsetDouble
Describe a linear system.
std::map< int, Composite > CompositeMap
Definition: MeshGraph.h:115
StdRegions::Orientation orient
Orientation of entity within higher dimensional entity.
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
void SetUpDG(const std::string="DefaultVar")
Set up all DG member variables and maps.
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
boost::shared_ptr< ElementEdgeVector > ElementEdgeVectorSharedPtr
Definition: MeshGraph.h:134
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:966
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
Definition: AssemblyMap.h:55
2D geometry information
Definition: Geometry2D.h:65
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:250
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1649
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
Definition: ExpList2D.h:60
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1493
virtual void v_GetFwdBwdTracePhys(const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
This method extracts the "forward" and "backward" trace data from the array field and puts the data i...
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
void EvaluateHDGPostProcessing(Array< OneD, NekDouble > &outarray)
Evaluate HDG post-processing to increase polynomial order of solution.
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
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Calculates the result of the multiplication of a global matrix of type specified by mkey with a vecto...
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:219
boost::shared_ptr< Basis > BasisSharedPtr
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
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
void GenerateBoundaryConditionExpansion(const SpatialDomains::MeshGraphSharedPtr &graph2D, const SpatialDomains::BoundaryConditions &bcs, const std::string &variable, const bool DeclareCoeffPhysArrays=true)
This function discretises the boundary conditions by setting up a list of one-dimensional boundary ex...
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Set up a list of element IDs and edge IDs that link to the boundary conditions.
#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
NekDouble L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error with respect to soln of the global spectral/hp element approximat...
Definition: ExpList.h:528
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
NekDouble L2_DGDeriv(const int dir, const Array< OneD, const NekDouble > &soln)
Calculate the error of the derivative using the consistent DG evaluation of .
Describes the specification for a Basis.
Definition: Basis.h:50
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:230
boost::shared_ptr< RobinBoundaryCondition > RobinBCShPtr
Definition: Conditions.h:222
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
boost::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:49