MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fespace.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 // Implementation of FiniteElementSpace
13 
14 #include "fem.hpp"
15 
16 #include <cmath>
17 #include <cstdarg>
18 #include <limits>
19 
20 using namespace std;
21 
22 namespace mfem
23 {
24 
25 template <> void Ordering::
26 DofsToVDofs<Ordering::byNODES>(int ndofs, int vdim, Array<int> &dofs)
27 {
28  // static method
29  int size = dofs.Size();
30  dofs.SetSize(size*vdim);
31  for (int vd = 1; vd < vdim; vd++)
32  {
33  for (int i = 0; i < size; i++)
34  {
35  dofs[i+size*vd] = Map<byNODES>(ndofs, vdim, dofs[i], vd);
36  }
37  }
38 }
39 
40 template <> void Ordering::
41 DofsToVDofs<Ordering::byVDIM>(int ndofs, int vdim, Array<int> &dofs)
42 {
43  // static method
44  int size = dofs.Size();
45  dofs.SetSize(size*vdim);
46  for (int vd = vdim-1; vd >= 0; vd--)
47  {
48  for (int i = 0; i < size; i++)
49  {
50  dofs[i+size*vd] = Map<byVDIM>(ndofs, vdim, dofs[i], vd);
51  }
52  }
53 }
54 
55 int FiniteElementSpace::GetOrder(int i) const
56 {
57  int GeomType = mesh->GetElementBaseGeometry(i);
58  return fec->FiniteElementForGeometry(GeomType)->GetOrder();
59 }
60 
61 int FiniteElementSpace::GetFaceOrder(int i) const
62 {
63  int GeomType = mesh->GetFaceBaseGeometry(i);
64  return fec->FiniteElementForGeometry(GeomType)->GetOrder();
65 }
66 
67 void FiniteElementSpace::DofsToVDofs (Array<int> &dofs, int ndofs) const
68 {
69  if (vdim == 1) { return; }
70  if (ndofs < 0) { ndofs = this->ndofs; }
71 
72  if (ordering == Ordering::byNODES)
73  {
74  Ordering::DofsToVDofs<Ordering::byNODES>(ndofs, vdim, dofs);
75  }
76  else
77  {
78  Ordering::DofsToVDofs<Ordering::byVDIM>(ndofs, vdim, dofs);
79  }
80 }
81 
82 void FiniteElementSpace::DofsToVDofs(int vd, Array<int> &dofs, int ndofs) const
83 {
84  if (vdim == 1) { return; }
85  if (ndofs < 0) { ndofs = this->ndofs; }
86 
87  if (ordering == Ordering::byNODES)
88  {
89  for (int i = 0; i < dofs.Size(); i++)
90  {
91  dofs[i] = Ordering::Map<Ordering::byNODES>(ndofs, vdim, dofs[i], vd);
92  }
93  }
94  else
95  {
96  for (int i = 0; i < dofs.Size(); i++)
97  {
98  dofs[i] = Ordering::Map<Ordering::byVDIM>(ndofs, vdim, dofs[i], vd);
99  }
100  }
101 }
102 
103 int FiniteElementSpace::DofToVDof(int dof, int vd, int ndofs) const
104 {
105  if (vdim == 1) { return dof; }
106  if (ndofs < 0) { ndofs = this->ndofs; }
107 
108  if (ordering == Ordering::byNODES)
109  {
110  return Ordering::Map<Ordering::byNODES>(ndofs, vdim, dof, vd);
111  }
112  else
113  {
114  return Ordering::Map<Ordering::byVDIM>(ndofs, vdim, dof, vd);
115  }
116 }
117 
118 // static function
119 void FiniteElementSpace::AdjustVDofs (Array<int> &vdofs)
120 {
121  int n = vdofs.Size(), *vdof = vdofs;
122  for (int i = 0; i < n; i++)
123  {
124  int j;
125  if ((j = vdof[i]) < 0)
126  {
127  vdof[i] = -1-j;
128  }
129  }
130 }
131 
132 void FiniteElementSpace::GetElementVDofs(int i, Array<int> &vdofs) const
133 {
134  GetElementDofs(i, vdofs);
135  DofsToVDofs(vdofs);
136 }
137 
138 void FiniteElementSpace::GetBdrElementVDofs(int i, Array<int> &vdofs) const
139 {
140  GetBdrElementDofs(i, vdofs);
141  DofsToVDofs(vdofs);
142 }
143 
144 void FiniteElementSpace::GetFaceVDofs(int i, Array<int> &vdofs) const
145 {
146  GetFaceDofs(i, vdofs);
147  DofsToVDofs(vdofs);
148 }
149 
150 void FiniteElementSpace::GetEdgeVDofs(int i, Array<int> &vdofs) const
151 {
152  GetEdgeDofs(i, vdofs);
153  DofsToVDofs(vdofs);
154 }
155 
156 void FiniteElementSpace::GetVertexVDofs(int i, Array<int> &vdofs) const
157 {
158  GetVertexDofs(i, vdofs);
159  DofsToVDofs(vdofs);
160 }
161 
162 void FiniteElementSpace::GetElementInteriorVDofs(int i, Array<int> &vdofs) const
163 {
164  GetElementInteriorDofs(i, vdofs);
165  DofsToVDofs(vdofs);
166 }
167 
168 void FiniteElementSpace::GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const
169 {
170  GetEdgeInteriorDofs(i, vdofs);
171  DofsToVDofs(vdofs);
172 }
173 
174 void FiniteElementSpace::BuildElementToDofTable() const
175 {
176  if (elem_dof) { return; }
177 
178  Table *el_dof = new Table;
179  Array<int> dofs;
180  el_dof -> MakeI (mesh -> GetNE());
181  for (int i = 0; i < mesh -> GetNE(); i++)
182  {
183  GetElementDofs (i, dofs);
184  el_dof -> AddColumnsInRow (i, dofs.Size());
185  }
186  el_dof -> MakeJ();
187  for (int i = 0; i < mesh -> GetNE(); i++)
188  {
189  GetElementDofs (i, dofs);
190  el_dof -> AddConnections (i, (int *)dofs, dofs.Size());
191  }
192  el_dof -> ShiftUpI();
193  elem_dof = el_dof;
194 }
195 
196 void FiniteElementSpace::RebuildElementToDofTable()
197 {
198  delete elem_dof;
199  elem_dof = NULL;
200  BuildElementToDofTable();
201 }
202 
203 void FiniteElementSpace::BuildDofToArrays()
204 {
205  if (dof_elem_array.Size()) { return; }
206 
207  BuildElementToDofTable();
208 
209  dof_elem_array.SetSize (ndofs);
210  dof_ldof_array.SetSize (ndofs);
211  dof_elem_array = -1;
212  for (int i = 0; i < mesh -> GetNE(); i++)
213  {
214  const int *dofs = elem_dof -> GetRow(i);
215  const int n = elem_dof -> RowSize(i);
216  for (int j = 0; j < n; j++)
217  {
218  if (dof_elem_array[dofs[j]] < 0)
219  {
220  dof_elem_array[dofs[j]] = i;
221  dof_ldof_array[dofs[j]] = j;
222  }
223  }
224  }
225 }
226 
227 static void mark_dofs(const Array<int> &dofs, Array<int> &mark_array)
228 {
229  for (int i = 0; i < dofs.Size(); i++)
230  {
231  int k = dofs[i];
232  if (k < 0) { k = -1 - k; }
233  mark_array[k] = -1;
234  }
235 }
236 
237 void FiniteElementSpace::GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
238  Array<int> &ess_vdofs) const
239 {
240  Array<int> vdofs;
241 
242  ess_vdofs.SetSize(GetVSize());
243  ess_vdofs = 0;
244 
245  for (int i = 0; i < GetNBE(); i++)
246  {
247  if (bdr_attr_is_ess[GetBdrAttribute(i)-1])
248  {
249  GetBdrElementVDofs(i, vdofs);
250  mark_dofs(vdofs, ess_vdofs);
251  }
252  }
253 
254  // mark possible hidden boundary edges in a non-conforming mesh, also
255  // local DOFs affected by boundary elements on other processors
256  if (mesh->ncmesh)
257  {
258  Array<int> bdr_verts, bdr_edges;
259  mesh->ncmesh->GetBoundaryClosure(bdr_attr_is_ess, bdr_verts, bdr_edges);
260 
261  for (int i = 0; i < bdr_verts.Size(); i++)
262  {
263  GetVertexVDofs(bdr_verts[i], vdofs);
264  mark_dofs(vdofs, ess_vdofs);
265  }
266  for (int i = 0; i < bdr_edges.Size(); i++)
267  {
268  GetEdgeVDofs(bdr_edges[i], vdofs);
269  mark_dofs(vdofs, ess_vdofs);
270  }
271  }
272 }
273 
274 void FiniteElementSpace::GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
275  Array<int> &ess_tdof_list)
276 {
277  Array<int> ess_vdofs, ess_tdofs;
278  GetEssentialVDofs(bdr_attr_is_ess, ess_vdofs);
279  const SparseMatrix *R = GetConformingRestriction();
280  if (!R)
281  {
282  ess_tdofs.MakeRef(ess_vdofs);
283  }
284  else
285  {
286  R->BooleanMult(ess_vdofs, ess_tdofs);
287  }
288  MarkerToList(ess_tdofs, ess_tdof_list);
289 }
290 
291 // static method
292 void FiniteElementSpace::MarkerToList(const Array<int> &marker,
293  Array<int> &list)
294 {
295  int num_marked = 0;
296  for (int i = 0; i < marker.Size(); i++)
297  {
298  if (marker[i]) { num_marked++; }
299  }
300  list.SetSize(0);
301  list.Reserve(num_marked);
302  for (int i = 0; i < marker.Size(); i++)
303  {
304  if (marker[i]) { list.Append(i); }
305  }
306 }
307 
308 // static method
309 void FiniteElementSpace::ListToMarker(const Array<int> &list, int marker_size,
310  Array<int> &marker, int mark_val)
311 {
312  marker.SetSize(marker_size);
313  marker = 0;
314  for (int i = 0; i < list.Size(); i++)
315  {
316  marker[list[i]] = mark_val;
317  }
318 }
319 
320 void FiniteElementSpace::ConvertToConformingVDofs(const Array<int> &dofs,
321  Array<int> &cdofs)
322 {
323  GetConformingProlongation();
324  if (cP) { cP->BooleanMultTranspose(dofs, cdofs); }
325  else { dofs.Copy(cdofs); }
326 }
327 
328 void FiniteElementSpace::ConvertFromConformingVDofs(const Array<int> &cdofs,
329  Array<int> &dofs)
330 {
331  GetConformingRestriction();
332  if (cR) { cR->BooleanMultTranspose(cdofs, dofs); }
333  else { cdofs.Copy(dofs); }
334 }
335 
336 SparseMatrix *
337 FiniteElementSpace::D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes)
338 {
339  int i, j;
340  Array<int> d_vdofs, c_vdofs;
341  SparseMatrix *R;
342 
343  R = new SparseMatrix (cfes -> GetVSize(), GetVSize());
344 
345  for (i = 0; i < mesh -> GetNE(); i++)
346  {
347  this -> GetElementVDofs (i, d_vdofs);
348  cfes -> GetElementVDofs (i, c_vdofs);
349 
350 #ifdef MFEM_DEBUG
351  if (d_vdofs.Size() != c_vdofs.Size())
352  {
353  mfem_error ("FiniteElementSpace::D2C_GlobalRestrictionMatrix (...)");
354  }
355 #endif
356 
357  for (j = 0; j < d_vdofs.Size(); j++)
358  {
359  R -> Set (c_vdofs[j], d_vdofs[j], 1.0);
360  }
361  }
362 
363  R -> Finalize();
364 
365  return R;
366 }
367 
368 SparseMatrix *
369 FiniteElementSpace::D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
370 {
371  int i, j;
372  Array<int> d_dofs, c_dofs;
373  SparseMatrix *R;
374 
375  R = new SparseMatrix (cfes -> GetNDofs(), ndofs);
376 
377  for (i = 0; i < mesh -> GetNE(); i++)
378  {
379  this -> GetElementDofs (i, d_dofs);
380  cfes -> GetElementDofs (i, c_dofs);
381 
382 #ifdef MFEM_DEBUG
383  if (c_dofs.Size() != 1)
384  mfem_error ("FiniteElementSpace::"
385  "D2Const_GlobalRestrictionMatrix (...)");
386 #endif
387 
388  for (j = 0; j < d_dofs.Size(); j++)
389  {
390  R -> Set (c_dofs[0], d_dofs[j], 1.0);
391  }
392  }
393 
394  R -> Finalize();
395 
396  return R;
397 }
398 
399 SparseMatrix *
400 FiniteElementSpace::H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes)
401 {
402  SparseMatrix *R;
403  DenseMatrix loc_restr;
404  Array<int> l_dofs, h_dofs;
405 
406  R = new SparseMatrix (lfes -> GetNDofs(), ndofs);
407 
408  if (!lfes->GetNE())
409  {
410  R->Finalize();
411  return R;
412  }
413 
414  const FiniteElement *h_fe = this -> GetFE (0);
415  const FiniteElement *l_fe = lfes -> GetFE (0);
418  h_fe->Project(*l_fe, T, loc_restr);
419 
420  for (int i = 0; i < mesh -> GetNE(); i++)
421  {
422  this -> GetElementDofs (i, h_dofs);
423  lfes -> GetElementDofs (i, l_dofs);
424 
425  R -> SetSubMatrix (l_dofs, h_dofs, loc_restr, 1);
426  }
427 
428  R -> Finalize();
429 
430  return R;
431 }
432 
433 static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
434  Array<int>& slave_dofs, DenseMatrix& I)
435 {
436  for (int i = 0; i < slave_dofs.Size(); i++)
437  {
438  int sdof = slave_dofs[i];
439  if (!deps.RowSize(sdof)) // not processed yet?
440  {
441  for (int j = 0; j < master_dofs.Size(); j++)
442  {
443  double coef = I(i, j);
444  if (std::abs(coef) > 1e-12)
445  {
446  int mdof = master_dofs[j];
447  if (mdof != sdof && mdof != (-1-sdof))
448  {
449  deps.Add(sdof, mdof, coef);
450  }
451  }
452  }
453  }
454  }
455 }
456 
457 static bool DofFinalizable(int dof, const Array<bool>& finalized,
458  const SparseMatrix& deps)
459 {
460  const int* dep = deps.GetRowColumns(dof);
461  int ndep = deps.RowSize(dof);
462 
463  // are all constraining DOFs finalized?
464  for (int i = 0; i < ndep; i++)
465  {
466  if (!finalized[dep[i]]) { return false; }
467  }
468  return true;
469 }
470 
474 void FiniteElementSpace::GetEdgeFaceDofs(int type, int index, Array<int> &dofs)
475 {
476  dofs.SetSize(0);
477  if (type)
478  {
479  if (index < mesh->GetNFaces()) { GetFaceDofs(index, dofs); }
480  }
481  else
482  {
483  if (index < mesh->GetNEdges()) { GetEdgeDofs(index, dofs); }
484  }
485 }
486 
487 void FiniteElementSpace::GetConformingInterpolation()
488 {
489 #ifdef MFEM_USE_MPI
490  MFEM_VERIFY(dynamic_cast<ParFiniteElementSpace*>(this) == NULL,
491  "This method should not be used with a ParFiniteElementSpace!");
492 #endif
493  if (cP_is_set) { return; }
494  cP_is_set = true;
495 
496  // For each slave DOF, the dependency matrix will contain a row that
497  // expresses the slave DOF as a linear combination of its immediate master
498  // DOFs. Rows of independent DOFs will remain empty.
499  SparseMatrix deps(ndofs);
500 
501  // collect local edge/face dependencies
502  for (int type = 0; type <= 1; type++)
503  {
504  const NCMesh::NCList &list = type ? mesh->ncmesh->GetFaceList()
505  : mesh->ncmesh->GetEdgeList();
506  if (!list.masters.size()) { continue; }
507 
509  if (type) { T.SetFE(&QuadrilateralFE); }
510  else { T.SetFE(&SegmentFE); }
511 
512  int geom = type ? Geometry::SQUARE : Geometry::SEGMENT;
513  const FiniteElement* fe = fec->FiniteElementForGeometry(geom);
514  if (!fe) { continue; }
515 
516  Array<int> master_dofs, slave_dofs;
517  DenseMatrix I(fe->GetDof());
518 
519  // loop through all master edges/faces, constrain their slave edges/faces
520  for (unsigned mi = 0; mi < list.masters.size(); mi++)
521  {
522  const NCMesh::Master &master = list.masters[mi];
523  GetEdgeFaceDofs(type, master.index, master_dofs);
524  if (!master_dofs.Size()) { continue; }
525 
526  for (int si = master.slaves_begin; si < master.slaves_end; si++)
527  {
528  const NCMesh::Slave &slave = list.slaves[si];
529  GetEdgeFaceDofs(type, slave.index, slave_dofs);
530  if (!slave_dofs.Size()) { continue; }
531 
532  slave.OrientedPointMatrix(T.GetPointMat());
533  fe->GetLocalInterpolation(T, I);
534 
535  // make each slave DOF dependent on all master DOFs
536  AddDependencies(deps, master_dofs, slave_dofs, I);
537  }
538  }
539  }
540 
541  deps.Finalize();
542 
543  // DOFs that stayed independent are true DOFs
544  int n_true_dofs = 0;
545  for (int i = 0; i < ndofs; i++)
546  {
547  if (!deps.RowSize(i)) { n_true_dofs++; }
548  }
549 
550  // if all dofs are true dofs leave cP and cR NULL
551  if (n_true_dofs == ndofs)
552  {
553  cP = cR = NULL; // will be treated as identities
554  return;
555  }
556 
557  // create the conforming restriction matrix cR
558  int *cR_J;
559  {
560  int *cR_I = new int[n_true_dofs+1];
561  double *cR_A = new double[n_true_dofs];
562  cR_J = new int[n_true_dofs];
563  for (int i = 0; i < n_true_dofs; i++)
564  {
565  cR_I[i] = i;
566  cR_A[i] = 1.0;
567  }
568  cR_I[n_true_dofs] = n_true_dofs;
569  cR = new SparseMatrix(cR_I, cR_J, cR_A, n_true_dofs, ndofs);
570  }
571 
572  // create the conforming prolongation matrix cP
573  cP = new SparseMatrix(ndofs, n_true_dofs);
574 
575  Array<bool> finalized(ndofs);
576  finalized = false;
577 
578  // put identity in the restriction and prolongation matrices for true DOFs
579  for (int i = 0, true_dof = 0; i < ndofs; i++)
580  {
581  if (!deps.RowSize(i))
582  {
583  cR_J[true_dof] = i;
584  cP->Add(i, true_dof++, 1.0);
585  finalized[i] = true;
586  }
587  }
588 
589  // Now calculate cP rows of slave DOFs as combinations of cP rows of their
590  // master DOFs. It is possible that some slave DOFs depend on DOFs that are
591  // themselves slaves. Here we resolve such indirect constraints by first
592  // calculating rows of the cP matrix for DOFs whose master DOF cP rows are
593  // already known (in the first iteration these are the true DOFs). In the
594  // second iteration, slaves of slaves can be 'finalized' (given a row in the
595  // cP matrix), in the third iteration slaves of slaves of slaves, etc.
596  bool finished;
597  int n_finalized = n_true_dofs;
598  Array<int> cols;
599  Vector srow;
600  do
601  {
602  finished = true;
603  for (int dof = 0; dof < ndofs; dof++)
604  {
605  if (!finalized[dof] && DofFinalizable(dof, finalized, deps))
606  {
607  const int* dep_col = deps.GetRowColumns(dof);
608  const double* dep_coef = deps.GetRowEntries(dof);
609  int n_dep = deps.RowSize(dof);
610 
611  for (int j = 0; j < n_dep; j++)
612  {
613  cP->GetRow(dep_col[j], cols, srow);
614  srow *= dep_coef[j];
615  cP->AddRow(dof, cols, srow);
616  }
617 
618  finalized[dof] = true;
619  n_finalized++;
620  finished = false;
621  }
622  }
623  }
624  while (!finished);
625 
626  // if everything is consistent (mesh, face orientations, etc.), we should
627  // be able to finalize all slave DOFs, otherwise it's a serious error
628  if (n_finalized != ndofs)
629  {
630  MFEM_ABORT("Error creating cP matrix.");
631  }
632 
633  cP->Finalize();
634 
635  if (vdim > 1)
636  {
637  MakeVDimMatrix(*cP);
638  MakeVDimMatrix(*cR);
639  }
640 }
641 
642 void FiniteElementSpace::MakeVDimMatrix(SparseMatrix &mat) const
643 {
644  if (vdim == 1) { return; }
645 
646  int height = mat.Height();
647  int width = mat.Width();
648 
649  SparseMatrix *vmat = new SparseMatrix(vdim*height, vdim*width);
650 
651  Array<int> dofs, vdofs;
652  Vector srow;
653  for (int i = 0; i < height; i++)
654  {
655  mat.GetRow(i, dofs, srow);
656  for (int vd = 0; vd < vdim; vd++)
657  {
658  dofs.Copy(vdofs);
659  DofsToVDofs(vd, vdofs, width);
660  vmat->SetRow(DofToVDof(i, vd, height), vdofs, srow);
661  }
662  }
663  vmat->Finalize();
664 
665  mat.Swap(*vmat);
666  delete vmat;
667 }
668 
669 const SparseMatrix* FiniteElementSpace::GetConformingProlongation()
670 {
671  if (Conforming()) { return NULL; }
672  if (!cP_is_set) { GetConformingInterpolation(); }
673  return cP;
674 }
675 
676 const SparseMatrix* FiniteElementSpace::GetConformingRestriction()
677 {
678  if (Conforming()) { return NULL; }
679  if (!cP_is_set) { GetConformingInterpolation(); }
680  return cR;
681 }
682 
683 int FiniteElementSpace::GetNConformingDofs()
684 {
685  const SparseMatrix* P = GetConformingProlongation();
686  return P ? (P->Width() / vdim) : ndofs;
687 }
688 
689 SparseMatrix* FiniteElementSpace::RefinementMatrix(int old_ndofs,
690  const Table* old_elem_dof)
691 {
692  MFEM_VERIFY(mesh->GetLastOperation() == Mesh::REFINE, "");
693  MFEM_VERIFY(ndofs >= old_ndofs, "Previous space is not coarser.");
694 
695  Array<int> dofs, old_dofs, old_vdofs;
696  Vector row;
697 
698  const CoarseFineTransformations &rtrans = mesh->GetRefinementTransforms();
699 
700  int geom = mesh->GetElementBaseGeometry(); // assuming the same geom
701  const FiniteElement *fe = fec->FiniteElementForGeometry(geom);
702 
704  isotr.SetIdentityTransformation(geom);
705 
706  int nmat = rtrans.point_matrices.SizeK();
707  int ldof = fe->GetDof(); // assuming the same FE everywhere
708 
709  // calculate local interpolation matrices for all refinement types
710  DenseTensor localP(ldof, ldof, nmat);
711  for (int i = 0; i < nmat; i++)
712  {
713  isotr.GetPointMat() = rtrans.point_matrices(i);
714  fe->GetLocalInterpolation(isotr, localP(i));
715  }
716 
717  SparseMatrix *P = new SparseMatrix(ndofs*vdim, old_ndofs*vdim, ldof);
718 
719  Array<int> mark(P->Height());
720  mark = 0;
721  for (int k = 0; k < mesh->GetNE(); k++)
722  {
723  const Embedding &emb = rtrans.embeddings[k];
724  DenseMatrix &lP = localP(emb.matrix);
725 
726  elem_dof->GetRow(k, dofs);
727  old_elem_dof->GetRow(emb.parent, old_dofs);
728 
729  for (int vd = 0; vd < vdim; vd++)
730  {
731  old_dofs.Copy(old_vdofs);
732  DofsToVDofs(vd, old_vdofs, old_ndofs);
733 
734  for (int i = 0; i < ldof; i++)
735  {
736  int r = DofToVDof(dofs[i], vd);
737  int m = (r >= 0) ? r : (-1 - r);
738 
739  if (!mark[m])
740  {
741  lP.GetRow(i, row);
742  P->SetRow(r, old_vdofs, row);
743  mark[m] = 1;
744  }
745  }
746  }
747  }
748 
749  MFEM_ASSERT(mark.Sum() == P->Height(), "Not all rows of P set.");
750  return P;
751 }
752 
754  const DenseMatrix &invdfdx,
755  const IntegrationPoint &pt, Vector &x)
756 {
757  // invert a linear transform with one Newton step
758  IntegrationPoint p0;
759  p0.Set3(0, 0, 0);
760  trans.Transform(p0, x);
761 
762  double store[3];
763  Vector v(store, x.Size());
764  pt.Get(v, x.Size());
765  v -= x;
766 
767  invdfdx.Mult(v, x);
768 }
769 
770 void FiniteElementSpace::GetLocalDerefinementMatrices(
771  int geom, const CoarseFineTransformations &dt, DenseTensor &localR)
772 {
773  const FiniteElement *fe = fec->FiniteElementForGeometry(geom);
774  const IntegrationRule &nodes = fe->GetNodes();
775 
776  LinearFECollection linfec;
778  isotr.SetFE(linfec.FiniteElementForGeometry(geom));
779 
780  int nmat = dt.point_matrices.SizeK();
781  int ldof = fe->GetDof();
782  int dim = mesh->Dimension();
783 
784  DenseMatrix invdfdx(dim);
785  IntegrationPoint ipt;
786  Vector pt(&ipt.x, dim), shape(ldof);
787 
788  // calculate local restriction matrices for all refinement types
789  localR.SetSize(ldof, ldof, nmat);
790  for (int i = 0; i < nmat; i++)
791  {
792  DenseMatrix &lR = localR(i);
793  lR = numeric_limits<double>::infinity(); // marks invalid rows
794 
795  isotr.GetPointMat() = dt.point_matrices(i);
796  isotr.SetIntPoint(&nodes[0]);
797  CalcInverse(isotr.Jacobian(), invdfdx);
798 
799  for (int j = 0; j < nodes.Size(); j++)
800  {
801  InvertLinearTrans(isotr, invdfdx, nodes[j], pt);
802  if (Geometries.CheckPoint(geom, ipt)) // do we need an epsilon here?
803  {
804  IntegrationPoint ip;
805  ip.Set(pt, dim);
806  fe->CalcShape(ip, shape); // TODO: H(curl), etc.?
807  MFEM_ASSERT(dynamic_cast<const NodalFiniteElement*>(fe),
808  "only nodal FEs are implemented");
809  lR.SetRow(j, shape);
810  }
811  }
812  }
813 }
814 
815 SparseMatrix* FiniteElementSpace::DerefinementMatrix(int old_ndofs,
816  const Table* old_elem_dof)
817 {
818  MFEM_VERIFY(Nonconforming(), "Not implemented for conforming meshes.");
819  MFEM_VERIFY(old_ndofs, "Missing previous (finer) space.");
820  MFEM_VERIFY(ndofs <= old_ndofs, "Previous space is not finer.");
821 
822  Array<int> dofs, old_dofs, old_vdofs;
823  Vector row;
824 
825  const CoarseFineTransformations &dtrans =
826  mesh->ncmesh->GetDerefinementTransforms();
827 
828  int geom = mesh->ncmesh->GetElementGeometry();
829  int ldof = fec->FiniteElementForGeometry(geom)->GetDof();
830 
831  DenseTensor localR;
832  GetLocalDerefinementMatrices(geom, dtrans, localR);
833 
834  SparseMatrix *R = new SparseMatrix(ndofs*vdim, old_ndofs*vdim, ldof);
835 
836  Array<int> mark(R->Height());
837  mark = 0;
838  for (int k = 0; k < dtrans.embeddings.Size(); k++)
839  {
840  const Embedding &emb = dtrans.embeddings[k];
841  DenseMatrix &lR = localR(emb.matrix);
842 
843  elem_dof->GetRow(emb.parent, dofs);
844  old_elem_dof->GetRow(k, old_dofs);
845 
846  for (int vd = 0; vd < vdim; vd++)
847  {
848  old_dofs.Copy(old_vdofs);
849  DofsToVDofs(vd, old_vdofs, old_ndofs);
850 
851  for (int i = 0; i < lR.Height(); i++)
852  {
853  if (lR(i, 0) == numeric_limits<double>::infinity()) { continue; }
854 
855  int r = DofToVDof(dofs[i], vd);
856  int m = (r >= 0) ? r : (-1 - r);
857 
858  if (!mark[m])
859  {
860  lR.GetRow(i, row);
861  R->SetRow(r, old_vdofs, row);
862  mark[m] = 1;
863  }
864  }
865  }
866  }
867 
868  MFEM_ASSERT(mark.Sum() == R->Height(), "Not all rows of R set.");
869  return R;
870 }
871 
872 FiniteElementSpace::FiniteElementSpace(Mesh *mesh,
873  const FiniteElementCollection *fec,
874  int vdim, int ordering)
875 {
876  this->mesh = mesh;
877  this->fec = fec;
878  this->vdim = vdim;
879  this->ordering = (Ordering::Type) ordering;
880 
881  elem_dof = NULL;
882  sequence = mesh->GetSequence();
883 
884  const NURBSFECollection *nurbs_fec =
885  dynamic_cast<const NURBSFECollection *>(fec);
886  if (nurbs_fec)
887  {
888  if (!mesh->NURBSext)
889  {
890  mfem_error("FiniteElementSpace::FiniteElementSpace :\n"
891  " NURBS FE space requires NURBS mesh.");
892  }
893  else
894  {
895  int Order = nurbs_fec->GetOrder();
896  if (mesh->NURBSext->GetOrder() == Order)
897  {
898  NURBSext = mesh->NURBSext;
899  own_ext = 0;
900  }
901  else
902  {
903  NURBSext = new NURBSExtension(mesh->NURBSext, Order);
904  own_ext = 1;
905  }
906  UpdateNURBS();
907  cP = cR = NULL;
908  cP_is_set = false;
909  T = NULL;
910  own_T = true;
911  }
912  }
913  else
914  {
915  NURBSext = NULL;
916  own_ext = 0;
917  Construct();
918  }
919 
920  BuildElementToDofTable();
921 }
922 
923 NURBSExtension *FiniteElementSpace::StealNURBSext()
924 {
925  if (NURBSext && !own_ext)
926  {
927  mfem_error("FiniteElementSpace::StealNURBSext");
928  }
929  own_ext = 0;
930 
931  return NURBSext;
932 }
933 
934 void FiniteElementSpace::UpdateNURBS()
935 {
936  nvdofs = 0;
937  nedofs = 0;
938  nfdofs = 0;
939  nbdofs = 0;
940  fdofs = NULL;
941  bdofs = NULL;
942 
943  dynamic_cast<const NURBSFECollection *>(fec)->Reset();
944 
945  ndofs = NURBSext->GetNDof();
946  elem_dof = NURBSext->GetElementDofTable();
947  bdrElem_dof = NURBSext->GetBdrElementDofTable();
948 }
949 
950 void FiniteElementSpace::Construct()
951 {
952  int i;
953 
954  elem_dof = NULL;
955  bdrElem_dof = NULL;
956 
957  nvdofs = mesh->GetNV() * fec->DofForGeometry(Geometry::POINT);
958 
959  if ( mesh->Dimension() > 1 )
960  {
961  nedofs = mesh->GetNEdges() * fec->DofForGeometry(Geometry::SEGMENT);
962  }
963  else
964  {
965  nedofs = 0;
966  }
967 
968  ndofs = 0;
969  nfdofs = 0;
970  nbdofs = 0;
971  bdofs = NULL;
972  fdofs = NULL;
973  cP = NULL;
974  cR = NULL;
975  cP_is_set = false;
976  T = NULL;
977  own_T = true;
978 
979  if (mesh->Dimension() == 3 && mesh->GetNE())
980  {
981  // Here we assume that all faces in the mesh have the same base
982  // geometry -- the base geometry of the 0-th face element.
983  // The class Mesh assumes the same inside GetFaceBaseGeometry(...).
984  // Thus we do not need to generate all the faces in the mesh
985  // if we do not need them.
986  int fdof = fec->DofForGeometry(mesh->GetFaceBaseGeometry(0));
987  if (fdof > 0)
988  {
989  fdofs = new int[mesh->GetNFaces()+1];
990  fdofs[0] = 0;
991  for (i = 0; i < mesh->GetNFaces(); i++)
992  {
993  nfdofs += fdof;
994  // nfdofs += fec->DofForGeometry(mesh->GetFaceBaseGeometry(i));
995  fdofs[i+1] = nfdofs;
996  }
997  }
998  }
999 
1000  bdofs = new int[mesh->GetNE()+1];
1001  bdofs[0] = 0;
1002  for (i = 0; i < mesh->GetNE(); i++)
1003  {
1004  nbdofs += fec->DofForGeometry(mesh->GetElementBaseGeometry(i));
1005  bdofs[i+1] = nbdofs;
1006  }
1007 
1008  ndofs = nvdofs + nedofs + nfdofs + nbdofs;
1009 
1010  // Do not build elem_dof Table here: in parallel it has to be constructed
1011  // later.
1012 }
1013 
1014 void FiniteElementSpace::GetElementDofs (int i, Array<int> &dofs) const
1015 {
1016  if (elem_dof)
1017  {
1018  elem_dof -> GetRow (i, dofs);
1019  }
1020  else
1021  {
1022  Array<int> V, E, Eo, F, Fo;
1023  int k, j, nv, ne, nf, nb, nfd, nd;
1024  int *ind, dim;
1025 
1026  dim = mesh->Dimension();
1027  nv = fec->DofForGeometry(Geometry::POINT);
1028  ne = (dim > 1) ? ( fec->DofForGeometry(Geometry::SEGMENT) ) : ( 0 );
1029  nb = fec->DofForGeometry(mesh->GetElementBaseGeometry(i));
1030  if (nv > 0)
1031  {
1032  mesh->GetElementVertices(i, V);
1033  }
1034  if (ne > 0)
1035  {
1036  mesh->GetElementEdges(i, E, Eo);
1037  }
1038  nfd = 0;
1039  if (dim == 3)
1040  {
1041  if (fec->HasFaceDofs(mesh->GetElementBaseGeometry(i)))
1042  {
1043  mesh->GetElementFaces(i, F, Fo);
1044  for (k = 0; k < F.Size(); k++)
1045  {
1046  nfd += fec->DofForGeometry(mesh->GetFaceBaseGeometry(F[k]));
1047  }
1048  }
1049  }
1050  nd = V.Size() * nv + E.Size() * ne + nfd + nb;
1051  dofs.SetSize(nd);
1052  if (nv > 0)
1053  {
1054  for (k = 0; k < V.Size(); k++)
1055  {
1056  for (j = 0; j < nv; j++)
1057  {
1058  dofs[k*nv+j] = V[k]*nv+j;
1059  }
1060  }
1061  nv *= V.Size();
1062  }
1063  if (ne > 0)
1064  {
1065  // if (dim > 1)
1066  for (k = 0; k < E.Size(); k++)
1067  {
1068  ind = fec->DofOrderForOrientation(Geometry::SEGMENT, Eo[k]);
1069  for (j = 0; j < ne; j++)
1070  {
1071  if (ind[j] < 0)
1072  {
1073  dofs[nv+k*ne+j] = -1 - ( nvdofs+E[k]*ne+(-1-ind[j]) );
1074  }
1075  else
1076  {
1077  dofs[nv+k*ne+j] = nvdofs+E[k]*ne+ind[j];
1078  }
1079  }
1080  }
1081  }
1082  ne = nv + ne * E.Size();
1083  if (nfd > 0)
1084  // if (dim == 3)
1085  {
1086  for (k = 0; k < F.Size(); k++)
1087  {
1088  ind = fec->DofOrderForOrientation(mesh->GetFaceBaseGeometry(F[k]),
1089  Fo[k]);
1090  nf = fec->DofForGeometry(mesh->GetFaceBaseGeometry(F[k]));
1091  for (j = 0; j < nf; j++)
1092  {
1093  if (ind[j] < 0)
1094  {
1095  dofs[ne+j] = -1 - ( nvdofs+nedofs+fdofs[F[k]]+(-1-ind[j]) );
1096  }
1097  else
1098  {
1099  dofs[ne+j] = nvdofs+nedofs+fdofs[F[k]]+ind[j];
1100  }
1101  }
1102  ne += nf;
1103  }
1104  }
1105  k = nvdofs + nedofs + nfdofs + bdofs[i];
1106  for (j = 0; j < nb; j++)
1107  {
1108  dofs[ne+j] = k + j;
1109  }
1110  }
1111 }
1112 
1113 const FiniteElement *FiniteElementSpace::GetFE(int i) const
1114 {
1115  const FiniteElement *FE =
1116  fec->FiniteElementForGeometry(mesh->GetElementBaseGeometry(i));
1117 
1118  if (NURBSext)
1119  {
1120  NURBSext->LoadFE(i, FE);
1121  }
1122 
1123  return FE;
1124 }
1125 
1126 void FiniteElementSpace::GetBdrElementDofs(int i, Array<int> &dofs) const
1127 {
1128  if (bdrElem_dof)
1129  {
1130  bdrElem_dof->GetRow(i, dofs);
1131  }
1132  else
1133  {
1134  Array<int> V, E, Eo;
1135  int k, j, nv, ne, nf, nd, iF, oF;
1136  int *ind, dim;
1137 
1138  dim = mesh->Dimension();
1139  nv = fec->DofForGeometry(Geometry::POINT);
1140  if (nv > 0)
1141  {
1142  mesh->GetBdrElementVertices(i, V);
1143  }
1144  ne = (dim > 1) ? ( fec->DofForGeometry(Geometry::SEGMENT) ) : ( 0 );
1145  if (ne > 0)
1146  {
1147  mesh->GetBdrElementEdges(i, E, Eo);
1148  }
1149  nd = V.Size() * nv + E.Size() * ne;
1150  nf = (dim == 3) ? (fec->DofForGeometry(
1151  mesh->GetBdrElementBaseGeometry(i))) : (0);
1152  if (nf > 0)
1153  {
1154  nd += nf;
1155  mesh->GetBdrElementFace(i, &iF, &oF);
1156  }
1157  dofs.SetSize(nd);
1158  if (nv > 0)
1159  {
1160  for (k = 0; k < V.Size(); k++)
1161  {
1162  for (j = 0; j < nv; j++)
1163  {
1164  dofs[k*nv+j] = V[k]*nv+j;
1165  }
1166  }
1167  nv *= V.Size();
1168  }
1169  if (ne > 0)
1170  {
1171  // if (dim > 1)
1172  for (k = 0; k < E.Size(); k++)
1173  {
1174  ind = fec->DofOrderForOrientation(Geometry::SEGMENT, Eo[k]);
1175  for (j = 0; j < ne; j++)
1176  {
1177  if (ind[j] < 0)
1178  {
1179  dofs[nv+k*ne+j] = -1 - ( nvdofs+E[k]*ne+(-1-ind[j]) );
1180  }
1181  else
1182  {
1183  dofs[nv+k*ne+j] = nvdofs+E[k]*ne+ind[j];
1184  }
1185  }
1186  }
1187  }
1188  if (nf > 0)
1189  // if (dim == 3)
1190  {
1191  ne = nv + ne * E.Size();
1192  ind = (fec->DofOrderForOrientation(
1193  mesh->GetBdrElementBaseGeometry(i), oF));
1194  for (j = 0; j < nf; j++)
1195  {
1196  if (ind[j] < 0)
1197  {
1198  dofs[ne+j] = -1 - ( nvdofs+nedofs+fdofs[iF]+(-1-ind[j]) );
1199  }
1200  else
1201  {
1202  dofs[ne+j] = nvdofs+nedofs+fdofs[iF]+ind[j];
1203  }
1204  }
1205  }
1206  }
1207 }
1208 
1209 void FiniteElementSpace::GetFaceDofs(int i, Array<int> &dofs) const
1210 {
1211  int j, k, nv, ne, nf, nd, dim = mesh->Dimension();
1212  Array<int> V, E, Eo;
1213  const int *ind;
1214 
1215  // for 1D, 2D and 3D faces
1216  nv = fec->DofForGeometry(Geometry::POINT);
1217  ne = (dim > 1) ? fec->DofForGeometry(Geometry::SEGMENT) : 0;
1218  if (nv > 0)
1219  {
1220  mesh->GetFaceVertices(i, V);
1221  }
1222  if (ne > 0)
1223  {
1224  mesh->GetFaceEdges(i, E, Eo);
1225  }
1226  nf = (fdofs) ? (fdofs[i+1]-fdofs[i]) : (0);
1227  nd = V.Size() * nv + E.Size() * ne + nf;
1228  dofs.SetSize(nd);
1229  if (nv > 0)
1230  {
1231  for (k = 0; k < V.Size(); k++)
1232  {
1233  for (j = 0; j < nv; j++)
1234  {
1235  dofs[k*nv+j] = V[k]*nv+j;
1236  }
1237  }
1238  }
1239  nv *= V.Size();
1240  if (ne > 0)
1241  {
1242  for (k = 0; k < E.Size(); k++)
1243  {
1244  ind = fec->DofOrderForOrientation(Geometry::SEGMENT, Eo[k]);
1245  for (j = 0; j < ne; j++)
1246  {
1247  if (ind[j] < 0)
1248  {
1249  dofs[nv+k*ne+j] = -1 - ( nvdofs+E[k]*ne+(-1-ind[j]) );
1250  }
1251  else
1252  {
1253  dofs[nv+k*ne+j] = nvdofs+E[k]*ne+ind[j];
1254  }
1255  }
1256  }
1257  }
1258  ne = nv + ne * E.Size();
1259  if (nf > 0)
1260  {
1261  for (j = nvdofs+nedofs+fdofs[i], k = 0; k < nf; j++, k++)
1262  {
1263  dofs[ne+k] = j;
1264  }
1265  }
1266 }
1267 
1268 void FiniteElementSpace::GetEdgeDofs(int i, Array<int> &dofs) const
1269 {
1270  int j, k, nv, ne;
1271  Array<int> V;
1272 
1273  nv = fec->DofForGeometry(Geometry::POINT);
1274  if (nv > 0)
1275  {
1276  mesh->GetEdgeVertices(i, V);
1277  }
1278  ne = fec->DofForGeometry(Geometry::SEGMENT);
1279  dofs.SetSize(2*nv+ne);
1280  if (nv > 0)
1281  {
1282  for (k = 0; k < 2; k++)
1283  {
1284  for (j = 0; j < nv; j++)
1285  {
1286  dofs[k*nv+j] = V[k]*nv+j;
1287  }
1288  }
1289  }
1290  nv *= 2;
1291  for (j = 0, k = nvdofs+i*ne; j < ne; j++, k++)
1292  {
1293  dofs[nv+j] = k;
1294  }
1295 }
1296 
1297 void FiniteElementSpace::GetVertexDofs(int i, Array<int> &dofs) const
1298 {
1299  int j, nv;
1300 
1301  nv = fec->DofForGeometry(Geometry::POINT);
1302  dofs.SetSize(nv);
1303  for (j = 0; j < nv; j++)
1304  {
1305  dofs[j] = i*nv+j;
1306  }
1307 }
1308 
1309 void FiniteElementSpace::GetElementInteriorDofs (int i, Array<int> &dofs) const
1310 {
1311  int j, k, nb;
1312  nb = fec -> DofForGeometry (mesh -> GetElementBaseGeometry (i));
1313  dofs.SetSize (nb);
1314  k = nvdofs + nedofs + nfdofs + bdofs[i];
1315  for (j = 0; j < nb; j++)
1316  {
1317  dofs[j] = k + j;
1318  }
1319 }
1320 
1321 void FiniteElementSpace::GetEdgeInteriorDofs (int i, Array<int> &dofs) const
1322 {
1323  int j, k, ne;
1324 
1325  ne = fec -> DofForGeometry (Geometry::SEGMENT);
1326  dofs.SetSize (ne);
1327  for (j = 0, k = nvdofs+i*ne; j < ne; j++, k++)
1328  {
1329  dofs[j] = k;
1330  }
1331 }
1332 
1333 const FiniteElement *FiniteElementSpace::GetBE (int i) const
1334 {
1335  const FiniteElement *BE;
1336 
1337  switch ( mesh->Dimension() )
1338  {
1339  case 1:
1340  BE = fec->FiniteElementForGeometry(Geometry::POINT);
1341  break;
1342  case 2:
1343  BE = fec->FiniteElementForGeometry(Geometry::SEGMENT);
1344  break;
1345  case 3:
1346  default:
1347  BE = fec->FiniteElementForGeometry(
1348  mesh->GetBdrElementBaseGeometry(i));
1349  }
1350 
1351  if (NURBSext)
1352  {
1353  NURBSext->LoadBE(i, BE);
1354  }
1355 
1356  return BE;
1357 }
1358 
1359 const FiniteElement *FiniteElementSpace::GetFaceElement(int i) const
1360 {
1361  const FiniteElement *fe;
1362 
1363  switch (mesh->Dimension())
1364  {
1365  case 1:
1366  fe = fec->FiniteElementForGeometry(Geometry::POINT);
1367  break;
1368  case 2:
1369  fe = fec->FiniteElementForGeometry(Geometry::SEGMENT);
1370  break;
1371  case 3:
1372  default:
1373  fe = fec->FiniteElementForGeometry(mesh->GetFaceBaseGeometry(i));
1374  }
1375 
1376  // if (NURBSext)
1377  // NURBSext->LoadFaceElement(i, fe);
1378 
1379  return fe;
1380 }
1381 
1382 const FiniteElement *FiniteElementSpace::GetEdgeElement(int i) const
1383 {
1384  return fec->FiniteElementForGeometry(Geometry::SEGMENT);
1385 }
1386 
1387 const FiniteElement *FiniteElementSpace::GetTraceElement(
1388  int i, int geom_type) const
1389 {
1390  return fec->TraceFiniteElementForGeometry(geom_type);
1391 }
1392 
1393 FiniteElementSpace::~FiniteElementSpace()
1394 {
1395  Destroy();
1396 }
1397 
1398 void FiniteElementSpace::Destroy()
1399 {
1400  delete cR;
1401  delete cP;
1402  if (own_T) { delete T; }
1403 
1404  dof_elem_array.DeleteAll();
1405  dof_ldof_array.DeleteAll();
1406 
1407  if (NURBSext)
1408  {
1409  if (own_ext) { delete NURBSext; }
1410  }
1411  else
1412  {
1413  delete elem_dof;
1414  delete bdrElem_dof;
1415 
1416  delete [] bdofs;
1417  delete [] fdofs;
1418  }
1419 }
1420 
1421 void FiniteElementSpace::Update(bool want_transform)
1422 {
1423  if (mesh->GetSequence() == sequence)
1424  {
1425  return; // mesh and space are in sync, no-op
1426  }
1427  if (want_transform && mesh->GetSequence() != sequence + 1)
1428  {
1429  MFEM_ABORT("Error in update sequence. Space needs to be updated after "
1430  "each mesh modification.");
1431  }
1432  sequence = mesh->GetSequence();
1433 
1434  if (NURBSext)
1435  {
1436  UpdateNURBS();
1437  return;
1438  }
1439 
1440  Table* old_elem_dof= NULL;
1441  int old_ndofs;
1442 
1443  // save old DOF table
1444  if (want_transform)
1445  {
1446  old_elem_dof = elem_dof;
1447  elem_dof = NULL;
1448  old_ndofs = ndofs;
1449  }
1450 
1451  Destroy();
1452  Construct(); // sets T to NULL, own_T to true
1453  BuildElementToDofTable();
1454 
1455  if (want_transform)
1456  {
1457  // calculate appropriate GridFunction transformation
1458  switch (mesh->GetLastOperation())
1459  {
1460  case Mesh::REFINE:
1461  {
1462  T = RefinementMatrix(old_ndofs, old_elem_dof);
1463  break;
1464  }
1465 
1466  case Mesh::DEREFINE:
1467  {
1468  GetConformingInterpolation();
1469  T = DerefinementMatrix(old_ndofs, old_elem_dof);
1470  if (cP && cR)
1471  {
1472  T = new TripleProductOperator(cP, cR, T, false, false, true);
1473  }
1474  break;
1475  }
1476 
1477  default:
1478  break; // T stays NULL
1479  }
1480  }
1481 
1482  delete old_elem_dof;
1483 }
1484 
1485 void FiniteElementSpace::Save(std::ostream &out) const
1486 {
1487  out << "FiniteElementSpace\n"
1488  << "FiniteElementCollection: " << fec->Name() << '\n'
1489  << "VDim: " << vdim << '\n'
1490  << "Ordering: " << ordering << '\n';
1491 }
1492 
1493 }
Abstract class for Finite Elements.
Definition: fe.hpp:44
int GetOrder() const
Definition: fe_coll.hpp:278
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:263
void Set(const double *p, const int dim)
Definition: intrules.hpp:32
int Size() const
Logical size of the array.
Definition: array.hpp:109
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:212
void Get(double *p, const int dim) const
Definition: intrules.hpp:45
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:1816
Class for integration rule.
Definition: intrules.hpp:83
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:75
void InvertLinearTrans(IsoparametricTransformation &trans, const DenseMatrix &invdfdx, const IntegrationPoint &pt, Vector &x)
Definition: fespace.cpp:753
const Geometry::Type geom
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:259
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, but treat all elements as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:595
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:167
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
int SizeK() const
Definition: densemat.hpp:595
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:35
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:160
void SetIdentityTransformation(int GeomType)
Definition: eltrans.cpp:29
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:306
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
int Size() const
Returns the size of the vector.
Definition: vector.hpp:86
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:179
void SetSize(int i, int j, int k)
Definition: densemat.hpp:597
void OrientedPointMatrix(DenseMatrix &oriented_matrix) const
Return the point matrix oriented according to the master and slave edges.
Definition: ncmesh.cpp:2059
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:69
T Sum()
Sum all entries.
Definition: array.cpp:133
virtual void Transform(const IntegrationPoint &, Vector &)
Definition: eltrans.cpp:129
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:273
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:182
Geometry Geometries
Definition: geom.cpp:544
long GetSequence() const
Definition: mesh.hpp:871
int dim
Definition: ex3.cpp:47
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Definition: sparsemat.cpp:1952
DenseTensor point_matrices
matrices for IsoparametricTransformation
Definition: ncmesh.hpp:53
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:376
std::vector< Master > masters
Definition: ncmesh.hpp:170
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
const IntegrationRule & GetNodes() const
Definition: fe.hpp:113
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition: array.hpp:122
Array< Embedding > embeddings
fine element positions in their parents
Definition: ncmesh.hpp:54
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:454
void SetRow(int r, const Vector &row)
Definition: densemat.cpp:2684
virtual void Finalize(int skip_zeros=1)
Definition: sparsemat.cpp:706
virtual const DenseMatrix & Jacobian()
Definition: eltrans.cpp:52
int slaves_end
slave faces
Definition: ncmesh.hpp:145
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:3012
int GetGeomType() const
Returns the geometry type:
Definition: fe.hpp:85
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:63
std::vector< Slave > slaves
Definition: ncmesh.hpp:171
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
void mfem_error(const char *msg)
Definition: error.cpp:23
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:331
void Update(Vector &x, int k, DenseMatrix &h, Vector &s, Array< Vector * > &v)
Definition: solvers.cpp:500
NURBSExtension * NURBSext
Definition: mesh.hpp:143
Class for integration point with weight.
Definition: intrules.hpp:25
General triple product operator x -&gt; A*B*C*x, with ownership of the factors.
Definition: operator.hpp:191
void SetFE(const FiniteElement *FE)
Definition: eltrans.hpp:81
void GetRow(int r, Vector &row)
Definition: densemat.cpp:2189
int parent
element index in the coarse mesh
Definition: ncmesh.hpp:44
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:1991
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:474
Vector data type.
Definition: vector.hpp:33
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:142
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Definition: geom.cpp:255
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:51
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:2990
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:2173
int index
Mesh number.
Definition: ncmesh.hpp:133
BiLinear2DFiniteElement QuadrilateralFE
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:569
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
Defines the position of a fine element within a coarse element.
Definition: ncmesh.hpp:42
int matrix
index into CoarseFineTransformations::point_matrices
Definition: ncmesh.hpp:45