MFEM  v4.0
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 "../general/text.hpp"
15 #include "../general/forall.hpp"
16 #include "../mesh/mesh_headers.hpp"
17 #include "fem.hpp"
18 
19 #include <cmath>
20 #include <cstdarg>
21 #include <limits>
22 
23 using namespace std;
24 
25 namespace mfem
26 {
27 
28 template <> void Ordering::
29 DofsToVDofs<Ordering::byNODES>(int ndofs, int vdim, Array<int> &dofs)
30 {
31  // static method
32  int size = dofs.Size();
33  dofs.SetSize(size*vdim);
34  for (int vd = 1; vd < vdim; vd++)
35  {
36  for (int i = 0; i < size; i++)
37  {
38  dofs[i+size*vd] = Map<byNODES>(ndofs, vdim, dofs[i], vd);
39  }
40  }
41 }
42 
43 template <> void Ordering::
44 DofsToVDofs<Ordering::byVDIM>(int ndofs, int vdim, Array<int> &dofs)
45 {
46  // static method
47  int size = dofs.Size();
48  dofs.SetSize(size*vdim);
49  for (int vd = vdim-1; vd >= 0; vd--)
50  {
51  for (int i = 0; i < size; i++)
52  {
53  dofs[i+size*vd] = Map<byVDIM>(ndofs, vdim, dofs[i], vd);
54  }
55  }
56 }
57 
58 
59 FiniteElementSpace::FiniteElementSpace()
60  : mesh(NULL), fec(NULL), vdim(0), ordering(Ordering::byNODES),
61  ndofs(0), nvdofs(0), nedofs(0), nfdofs(0), nbdofs(0),
62  fdofs(NULL), bdofs(NULL),
63  elem_dof(NULL), bdrElem_dof(NULL),
64  NURBSext(NULL), own_ext(false),
65  cP(NULL), cR(NULL), cP_is_set(false),
66  Th(Operator::ANY_TYPE),
67  sequence(0)
68 { }
69 
71  Mesh *mesh,
72  const FiniteElementCollection *fec)
73 {
74  mesh = mesh ? mesh : orig.mesh;
75  fec = fec ? fec : orig.fec;
76  NURBSExtension *NURBSext = NULL;
77  if (orig.NURBSext && orig.NURBSext != orig.mesh->NURBSext)
78  {
79 #ifdef MFEM_USE_MPI
80  ParNURBSExtension *pNURBSext =
81  dynamic_cast<ParNURBSExtension *>(orig.NURBSext);
82  if (pNURBSext)
83  {
84  NURBSext = new ParNURBSExtension(*pNURBSext);
85  }
86  else
87 #endif
88  {
89  NURBSext = new NURBSExtension(*orig.NURBSext);
90  }
91  }
92  Constructor(mesh, NURBSext, fec, orig.vdim, orig.ordering);
93 }
94 
96 {
98  return fec->FiniteElementForGeometry(GeomType)->GetOrder();
99 }
100 
102 {
103  Geometry::Type GeomType = mesh->GetFaceBaseGeometry(i);
104  return fec->FiniteElementForGeometry(GeomType)->GetOrder();
105 }
106 
107 void FiniteElementSpace::DofsToVDofs (Array<int> &dofs, int ndofs) const
108 {
109  if (vdim == 1) { return; }
110  if (ndofs < 0) { ndofs = this->ndofs; }
111 
113  {
114  Ordering::DofsToVDofs<Ordering::byNODES>(ndofs, vdim, dofs);
115  }
116  else
117  {
118  Ordering::DofsToVDofs<Ordering::byVDIM>(ndofs, vdim, dofs);
119  }
120 }
121 
122 void FiniteElementSpace::DofsToVDofs(int vd, Array<int> &dofs, int ndofs) const
123 {
124  if (vdim == 1) { return; }
125  if (ndofs < 0) { ndofs = this->ndofs; }
126 
128  {
129  for (int i = 0; i < dofs.Size(); i++)
130  {
131  dofs[i] = Ordering::Map<Ordering::byNODES>(ndofs, vdim, dofs[i], vd);
132  }
133  }
134  else
135  {
136  for (int i = 0; i < dofs.Size(); i++)
137  {
138  dofs[i] = Ordering::Map<Ordering::byVDIM>(ndofs, vdim, dofs[i], vd);
139  }
140  }
141 }
142 
143 int FiniteElementSpace::DofToVDof(int dof, int vd, int ndofs) const
144 {
145  if (vdim == 1) { return dof; }
146  if (ndofs < 0) { ndofs = this->ndofs; }
147 
149  {
150  return Ordering::Map<Ordering::byNODES>(ndofs, vdim, dof, vd);
151  }
152  else
153  {
154  return Ordering::Map<Ordering::byVDIM>(ndofs, vdim, dof, vd);
155  }
156 }
157 
158 // static function
160 {
161  int n = vdofs.Size(), *vdof = vdofs;
162  for (int i = 0; i < n; i++)
163  {
164  int j;
165  if ((j = vdof[i]) < 0)
166  {
167  vdof[i] = -1-j;
168  }
169  }
170 }
171 
173 {
174  GetElementDofs(i, vdofs);
175  DofsToVDofs(vdofs);
176 }
177 
179 {
180  GetBdrElementDofs(i, vdofs);
181  DofsToVDofs(vdofs);
182 }
183 
185 {
186  GetFaceDofs(i, vdofs);
187  DofsToVDofs(vdofs);
188 }
189 
191 {
192  GetEdgeDofs(i, vdofs);
193  DofsToVDofs(vdofs);
194 }
195 
197 {
198  GetVertexDofs(i, vdofs);
199  DofsToVDofs(vdofs);
200 }
201 
203 {
204  GetElementInteriorDofs(i, vdofs);
205  DofsToVDofs(vdofs);
206 }
207 
209 {
210  GetEdgeInteriorDofs(i, vdofs);
211  DofsToVDofs(vdofs);
212 }
213 
215 {
216  if (elem_dof) { return; }
217 
218  Table *el_dof = new Table;
219  Array<int> dofs;
220  el_dof -> MakeI (mesh -> GetNE());
221  for (int i = 0; i < mesh -> GetNE(); i++)
222  {
223  GetElementDofs (i, dofs);
224  el_dof -> AddColumnsInRow (i, dofs.Size());
225  }
226  el_dof -> MakeJ();
227  for (int i = 0; i < mesh -> GetNE(); i++)
228  {
229  GetElementDofs (i, dofs);
230  el_dof -> AddConnections (i, (int *)dofs, dofs.Size());
231  }
232  el_dof -> ShiftUpI();
233  elem_dof = el_dof;
234 }
235 
237 {
238  delete elem_dof;
239  elem_dof = NULL;
241 }
242 
244 {
245  Array<int> dof_marker(ndofs);
246 
247  dof_marker = -1;
248 
249  int *J = elem_dof->GetJ(), nnz = elem_dof->Size_of_connections();
250  for (int k = 0, dof_counter = 0; k < nnz; k++)
251  {
252  const int sdof = J[k]; // signed dof
253  const int dof = (sdof < 0) ? -1-sdof : sdof;
254  int new_dof = dof_marker[dof];
255  if (new_dof < 0)
256  {
257  dof_marker[dof] = new_dof = dof_counter++;
258  }
259  J[k] = (sdof < 0) ? -1-new_dof : new_dof; // preserve the sign of sdof
260  }
261 }
262 
264 {
265  if (dof_elem_array.Size()) { return; }
266 
268 
271  dof_elem_array = -1;
272  for (int i = 0; i < mesh -> GetNE(); i++)
273  {
274  const int *dofs = elem_dof -> GetRow(i);
275  const int n = elem_dof -> RowSize(i);
276  for (int j = 0; j < n; j++)
277  {
278  if (dof_elem_array[dofs[j]] < 0)
279  {
280  dof_elem_array[dofs[j]] = i;
281  dof_ldof_array[dofs[j]] = j;
282  }
283  }
284  }
285 }
286 
287 static void mark_dofs(const Array<int> &dofs, Array<int> &mark_array)
288 {
289  for (int i = 0; i < dofs.Size(); i++)
290  {
291  int k = dofs[i];
292  if (k < 0) { k = -1 - k; }
293  mark_array[k] = -1;
294  }
295 }
296 
298  Array<int> &ess_vdofs,
299  int component) const
300 {
301  Array<int> vdofs, dofs;
302 
303  ess_vdofs.SetSize(GetVSize());
304  ess_vdofs = 0;
305 
306  for (int i = 0; i < GetNBE(); i++)
307  {
308  if (bdr_attr_is_ess[GetBdrAttribute(i)-1])
309  {
310  if (component < 0)
311  {
312  // Mark all components.
313  GetBdrElementVDofs(i, vdofs);
314  mark_dofs(vdofs, ess_vdofs);
315  }
316  else
317  {
318  GetBdrElementDofs(i, dofs);
319  for (int d = 0; d < dofs.Size(); d++)
320  { dofs[d] = DofToVDof(dofs[d], component); }
321  mark_dofs(dofs, ess_vdofs);
322  }
323  }
324  }
325 
326  // mark possible hidden boundary edges in a non-conforming mesh, also
327  // local DOFs affected by boundary elements on other processors
328  if (mesh->ncmesh)
329  {
330  Array<int> bdr_verts, bdr_edges;
331  mesh->ncmesh->GetBoundaryClosure(bdr_attr_is_ess, bdr_verts, bdr_edges);
332 
333  for (int i = 0; i < bdr_verts.Size(); i++)
334  {
335  if (component < 0)
336  {
337  GetVertexVDofs(bdr_verts[i], vdofs);
338  mark_dofs(vdofs, ess_vdofs);
339  }
340  else
341  {
342  GetVertexDofs(bdr_verts[i], dofs);
343  for (int d = 0; d < dofs.Size(); d++)
344  { dofs[d] = DofToVDof(dofs[d], component); }
345  mark_dofs(dofs, ess_vdofs);
346  }
347  }
348  for (int i = 0; i < bdr_edges.Size(); i++)
349  {
350  if (component < 0)
351  {
352  GetEdgeVDofs(bdr_edges[i], vdofs);
353  mark_dofs(vdofs, ess_vdofs);
354  }
355  else
356  {
357  GetEdgeDofs(bdr_edges[i], dofs);
358  for (int d = 0; d < dofs.Size(); d++)
359  { dofs[d] = DofToVDof(dofs[d], component); }
360  mark_dofs(dofs, ess_vdofs);
361  }
362  }
363  }
364 }
365 
367  Array<int> &ess_tdof_list,
368  int component)
369 {
370  Array<int> ess_vdofs, ess_tdofs;
371  GetEssentialVDofs(bdr_attr_is_ess, ess_vdofs, component);
373  if (!R)
374  {
375  ess_tdofs.MakeRef(ess_vdofs);
376  }
377  else
378  {
379  R->BooleanMult(ess_vdofs, ess_tdofs);
380  }
381  MarkerToList(ess_tdofs, ess_tdof_list);
382 }
383 
384 // static method
386  Array<int> &list)
387 {
388  int num_marked = 0;
389  marker.HostRead(); // make sure we can read the array on host
390  for (int i = 0; i < marker.Size(); i++)
391  {
392  if (marker[i]) { num_marked++; }
393  }
394  list.SetSize(0);
395  list.Reserve(num_marked);
396  for (int i = 0; i < marker.Size(); i++)
397  {
398  if (marker[i]) { list.Append(i); }
399  }
400 }
401 
402 // static method
403 void FiniteElementSpace::ListToMarker(const Array<int> &list, int marker_size,
404  Array<int> &marker, int mark_val)
405 {
406  marker.SetSize(marker_size);
407  marker = 0;
408  for (int i = 0; i < list.Size(); i++)
409  {
410  marker[list[i]] = mark_val;
411  }
412 }
413 
415  Array<int> &cdofs)
416 {
418  if (cP) { cP->BooleanMultTranspose(dofs, cdofs); }
419  else { dofs.Copy(cdofs); }
420 }
421 
423  Array<int> &dofs)
424 {
426  if (cR) { cR->BooleanMultTranspose(cdofs, dofs); }
427  else { cdofs.Copy(dofs); }
428 }
429 
430 SparseMatrix *
432 {
433  int i, j;
434  Array<int> d_vdofs, c_vdofs;
435  SparseMatrix *R;
436 
437  R = new SparseMatrix (cfes -> GetVSize(), GetVSize());
438 
439  for (i = 0; i < mesh -> GetNE(); i++)
440  {
441  this -> GetElementVDofs (i, d_vdofs);
442  cfes -> GetElementVDofs (i, c_vdofs);
443 
444 #ifdef MFEM_DEBUG
445  if (d_vdofs.Size() != c_vdofs.Size())
446  {
447  mfem_error ("FiniteElementSpace::D2C_GlobalRestrictionMatrix (...)");
448  }
449 #endif
450 
451  for (j = 0; j < d_vdofs.Size(); j++)
452  {
453  R -> Set (c_vdofs[j], d_vdofs[j], 1.0);
454  }
455  }
456 
457  R -> Finalize();
458 
459  return R;
460 }
461 
462 SparseMatrix *
464 {
465  int i, j;
466  Array<int> d_dofs, c_dofs;
467  SparseMatrix *R;
468 
469  R = new SparseMatrix (cfes -> GetNDofs(), ndofs);
470 
471  for (i = 0; i < mesh -> GetNE(); i++)
472  {
473  this -> GetElementDofs (i, d_dofs);
474  cfes -> GetElementDofs (i, c_dofs);
475 
476 #ifdef MFEM_DEBUG
477  if (c_dofs.Size() != 1)
478  mfem_error ("FiniteElementSpace::"
479  "D2Const_GlobalRestrictionMatrix (...)");
480 #endif
481 
482  for (j = 0; j < d_dofs.Size(); j++)
483  {
484  R -> Set (c_dofs[0], d_dofs[j], 1.0);
485  }
486  }
487 
488  R -> Finalize();
489 
490  return R;
491 }
492 
493 SparseMatrix *
495 {
496  SparseMatrix *R;
497  DenseMatrix loc_restr;
498  Array<int> l_dofs, h_dofs;
499 
500  R = new SparseMatrix (lfes -> GetNDofs(), ndofs);
501 
502  Geometry::Type cached_geom = Geometry::INVALID;
503  const FiniteElement *h_fe = NULL;
504  const FiniteElement *l_fe = NULL;
506 
507  for (int i = 0; i < mesh -> GetNE(); i++)
508  {
509  this -> GetElementDofs (i, h_dofs);
510  lfes -> GetElementDofs (i, l_dofs);
511 
512  // Assuming 'loc_restr' depends only on the Geometry::Type.
514  if (geom != cached_geom)
515  {
516  h_fe = this -> GetFE (i);
517  l_fe = lfes -> GetFE (i);
519  h_fe->Project(*l_fe, T, loc_restr);
520  cached_geom = geom;
521  }
522 
523  R -> SetSubMatrix (l_dofs, h_dofs, loc_restr, 1);
524  }
525 
526  R -> Finalize();
527 
528  return R;
529 }
530 
531 void
533  Array<int>& slave_dofs, DenseMatrix& I)
534 {
535  for (int i = 0; i < slave_dofs.Size(); i++)
536  {
537  int sdof = slave_dofs[i];
538  if (!deps.RowSize(sdof)) // not processed yet?
539  {
540  for (int j = 0; j < master_dofs.Size(); j++)
541  {
542  double coef = I(i, j);
543  if (std::abs(coef) > 1e-12)
544  {
545  int mdof = master_dofs[j];
546  if (mdof != sdof && mdof != (-1-sdof))
547  {
548  deps.Add(sdof, mdof, coef);
549  }
550  }
551  }
552  }
553  }
554 }
555 
556 bool FiniteElementSpace::DofFinalizable(int dof, const Array<bool>& finalized,
557  const SparseMatrix& deps)
558 {
559  const int* dep = deps.GetRowColumns(dof);
560  int ndep = deps.RowSize(dof);
561 
562  // are all constraining DOFs finalized?
563  for (int i = 0; i < ndep; i++)
564  {
565  if (!finalized[dep[i]]) { return false; }
566  }
567  return true;
568 }
569 
570 void
571 FiniteElementSpace::GetEntityDofs(int entity, int index, Array<int> &dofs) const
572 {
573  switch (entity)
574  {
575  case 0: GetVertexDofs(index, dofs); break;
576  case 1: GetEdgeDofs(index, dofs); break;
577  case 2: GetFaceDofs(index, dofs); break;
578  }
579 }
580 
581 
583 {
584 #ifdef MFEM_USE_MPI
585  MFEM_VERIFY(dynamic_cast<const ParFiniteElementSpace*>(this) == NULL,
586  "This method should not be used with a ParFiniteElementSpace!");
587 #endif
588 
589  if (cP_is_set) { return; }
590  cP_is_set = true;
591 
592  // For each slave DOF, the dependency matrix will contain a row that
593  // expresses the slave DOF as a linear combination of its immediate master
594  // DOFs. Rows of independent DOFs will remain empty.
595  SparseMatrix deps(ndofs);
596 
597  // collect local edge/face dependencies
598  for (int entity = 1; entity <= 2; entity++)
599  {
600  const NCMesh::NCList &list = (entity > 1) ? mesh->ncmesh->GetFaceList()
601  /* */ : mesh->ncmesh->GetEdgeList();
602  if (!list.masters.size()) { continue; }
603 
605  if (entity > 1) { T.SetFE(&QuadrilateralFE); }
606  else { T.SetFE(&SegmentFE); }
607 
609  const FiniteElement* fe = fec->FiniteElementForGeometry(geom);
610  if (!fe) { continue; }
611 
612  Array<int> master_dofs, slave_dofs;
613  DenseMatrix I(fe->GetDof());
614 
615  // loop through all master edges/faces, constrain their slave edges/faces
616  for (unsigned mi = 0; mi < list.masters.size(); mi++)
617  {
618  const NCMesh::Master &master = list.masters[mi];
619  GetEntityDofs(entity, master.index, master_dofs);
620  if (!master_dofs.Size()) { continue; }
621 
622  for (int si = master.slaves_begin; si < master.slaves_end; si++)
623  {
624  const NCMesh::Slave &slave = list.slaves[si];
625  GetEntityDofs(entity, slave.index, slave_dofs);
626  if (!slave_dofs.Size()) { continue; }
627 
628  slave.OrientedPointMatrix(T.GetPointMat());
630  fe->GetLocalInterpolation(T, I);
631 
632  // make each slave DOF dependent on all master DOFs
633  AddDependencies(deps, master_dofs, slave_dofs, I);
634  }
635  }
636  }
637 
638  deps.Finalize();
639 
640  // DOFs that stayed independent are true DOFs
641  int n_true_dofs = 0;
642  for (int i = 0; i < ndofs; i++)
643  {
644  if (!deps.RowSize(i)) { n_true_dofs++; }
645  }
646 
647  // if all dofs are true dofs leave cP and cR NULL
648  if (n_true_dofs == ndofs)
649  {
650  cP = cR = NULL; // will be treated as identities
651  return;
652  }
653 
654  // create the conforming restriction matrix cR
655  int *cR_J;
656  {
657  int *cR_I = new int[n_true_dofs+1];
658  double *cR_A = new double[n_true_dofs];
659  cR_J = new int[n_true_dofs];
660  for (int i = 0; i < n_true_dofs; i++)
661  {
662  cR_I[i] = i;
663  cR_A[i] = 1.0;
664  }
665  cR_I[n_true_dofs] = n_true_dofs;
666  cR = new SparseMatrix(cR_I, cR_J, cR_A, n_true_dofs, ndofs);
667  }
668 
669  // create the conforming prolongation matrix cP
670  cP = new SparseMatrix(ndofs, n_true_dofs);
671 
672  Array<bool> finalized(ndofs);
673  finalized = false;
674 
675  // put identity in the restriction and prolongation matrices for true DOFs
676  for (int i = 0, true_dof = 0; i < ndofs; i++)
677  {
678  if (!deps.RowSize(i))
679  {
680  cR_J[true_dof] = i;
681  cP->Add(i, true_dof++, 1.0);
682  finalized[i] = true;
683  }
684  }
685 
686  // Now calculate cP rows of slave DOFs as combinations of cP rows of their
687  // master DOFs. It is possible that some slave DOFs depend on DOFs that are
688  // themselves slaves. Here we resolve such indirect constraints by first
689  // calculating rows of the cP matrix for DOFs whose master DOF cP rows are
690  // already known (in the first iteration these are the true DOFs). In the
691  // second iteration, slaves of slaves can be 'finalized' (given a row in the
692  // cP matrix), in the third iteration slaves of slaves of slaves, etc.
693  bool finished;
694  int n_finalized = n_true_dofs;
695  Array<int> cols;
696  Vector srow;
697  do
698  {
699  finished = true;
700  for (int dof = 0; dof < ndofs; dof++)
701  {
702  if (!finalized[dof] && DofFinalizable(dof, finalized, deps))
703  {
704  const int* dep_col = deps.GetRowColumns(dof);
705  const double* dep_coef = deps.GetRowEntries(dof);
706  int n_dep = deps.RowSize(dof);
707 
708  for (int j = 0; j < n_dep; j++)
709  {
710  cP->GetRow(dep_col[j], cols, srow);
711  srow *= dep_coef[j];
712  cP->AddRow(dof, cols, srow);
713  }
714 
715  finalized[dof] = true;
716  n_finalized++;
717  finished = false;
718  }
719  }
720  }
721  while (!finished);
722 
723  // if everything is consistent (mesh, face orientations, etc.), we should
724  // be able to finalize all slave DOFs, otherwise it's a serious error
725  if (n_finalized != ndofs)
726  {
727  MFEM_ABORT("Error creating cP matrix.");
728  }
729 
730  cP->Finalize();
731 
732  if (vdim > 1)
733  {
734  MakeVDimMatrix(*cP);
735  MakeVDimMatrix(*cR);
736  }
737 
738  if (Device::IsEnabled()) { cP->BuildTranspose(); }
739 }
740 
742 {
743  if (vdim == 1) { return; }
744 
745  int height = mat.Height();
746  int width = mat.Width();
747 
748  SparseMatrix *vmat = new SparseMatrix(vdim*height, vdim*width);
749 
750  Array<int> dofs, vdofs;
751  Vector srow;
752  for (int i = 0; i < height; i++)
753  {
754  mat.GetRow(i, dofs, srow);
755  for (int vd = 0; vd < vdim; vd++)
756  {
757  dofs.Copy(vdofs);
758  DofsToVDofs(vd, vdofs, width);
759  vmat->SetRow(DofToVDof(i, vd, height), vdofs, srow);
760  }
761  }
762  vmat->Finalize();
763 
764  mat.Swap(*vmat);
765  delete vmat;
766 }
767 
768 
770 {
771  if (Conforming()) { return NULL; }
773  return cP;
774 }
775 
777 {
778  if (Conforming()) { return NULL; }
780  return cR;
781 }
782 
784 {
786  return P ? (P->Width() / vdim) : ndofs;
787 }
788 
790  ElementDofOrdering e_ordering) const
791 {
792  // Check if we have a discontinuous space using the FE collection:
793  const L2_FECollection *dg_space = dynamic_cast<const L2_FECollection*>(fec);
794  if (dg_space) { return NULL; }
795  // TODO: support other DG collections.
796  if (e_ordering == ElementDofOrdering::LEXICOGRAPHIC)
797  {
798  if (L2E_lex.Ptr() == NULL)
799  {
800  L2E_lex.Reset(new ElementRestriction(*this, e_ordering));
801  }
802  return L2E_lex.Ptr();
803  }
804  // e_ordering == ElementDofOrdering::NATIVE
805  if (L2E_nat.Ptr() == NULL)
806  {
807  L2E_nat.Reset(new ElementRestriction(*this, e_ordering));
808  }
809  return L2E_nat.Ptr();
810 }
811 
813  const IntegrationRule &ir) const
814 {
815  for (int i = 0; i < E2Q_array.Size(); i++)
816  {
817  const QuadratureInterpolator *qi = E2Q_array[i];
818  if (qi->IntRule == &ir) { return qi; }
819  }
820 
821  QuadratureInterpolator *qi = new QuadratureInterpolator(*this, ir);
822  E2Q_array.Append(qi);
823  return qi;
824 }
825 
827  const QuadratureSpace &qs) const
828 {
829  for (int i = 0; i < E2Q_array.Size(); i++)
830  {
831  const QuadratureInterpolator *qi = E2Q_array[i];
832  if (qi->qspace == &qs) { return qi; }
833  }
834 
835  QuadratureInterpolator *qi = new QuadratureInterpolator(*this, qs);
836  E2Q_array.Append(qi);
837  return qi;
838 }
839 
841  const int coarse_ndofs, const Table &coarse_elem_dof,
842  const DenseTensor localP[]) const
843 {
844  MFEM_VERIFY(mesh->GetLastOperation() == Mesh::REFINE, "");
845 
846  Array<int> dofs, coarse_dofs, coarse_vdofs;
847  Vector row;
848 
849  Mesh::GeometryList elem_geoms(*mesh);
850 
851  SparseMatrix *P;
852  if (elem_geoms.Size() == 1)
853  {
854  const int coarse_ldof = localP[elem_geoms[0]].SizeJ();
855  P = new SparseMatrix(GetVSize(), coarse_ndofs*vdim, coarse_ldof);
856  }
857  else
858  {
859  P = new SparseMatrix(GetVSize(), coarse_ndofs*vdim);
860  }
861 
862  Array<int> mark(P->Height());
863  mark = 0;
864 
866 
867  for (int k = 0; k < mesh->GetNE(); k++)
868  {
869  const Embedding &emb = rtrans.embeddings[k];
871  const DenseMatrix &lP = localP[geom](emb.matrix);
872  const int fine_ldof = localP[geom].SizeI();
873 
874  elem_dof->GetRow(k, dofs);
875  coarse_elem_dof.GetRow(emb.parent, coarse_dofs);
876 
877  for (int vd = 0; vd < vdim; vd++)
878  {
879  coarse_dofs.Copy(coarse_vdofs);
880  DofsToVDofs(vd, coarse_vdofs, coarse_ndofs);
881 
882  for (int i = 0; i < fine_ldof; i++)
883  {
884  int r = DofToVDof(dofs[i], vd);
885  int m = (r >= 0) ? r : (-1 - r);
886 
887  if (!mark[m])
888  {
889  lP.GetRow(i, row);
890  P->SetRow(r, coarse_vdofs, row);
891  mark[m] = 1;
892  }
893  }
894  }
895  }
896 
897  MFEM_ASSERT(mark.Sum() == P->Height(), "Not all rows of P set.");
898  if (elem_geoms.Size() != 1) { P->Finalize(); }
899  return P;
900 }
901 
903  Geometry::Type geom, DenseTensor &localP) const
904 {
905  const FiniteElement *fe = fec->FiniteElementForGeometry(geom);
906 
908  const DenseTensor &pmats = rtrans.GetPointMatrices(geom);
909 
910  int nmat = pmats.SizeK();
911  int ldof = fe->GetDof(); // assuming the same FE everywhere
912 
914  isotr.SetIdentityTransformation(geom);
915 
916  // calculate local interpolation matrices for all refinement types
917  localP.SetSize(ldof, ldof, nmat);
918  for (int i = 0; i < nmat; i++)
919  {
920  isotr.GetPointMat() = pmats(i);
921  isotr.FinalizeTransformation();
922  fe->GetLocalInterpolation(isotr, localP(i));
923  }
924 }
925 
927  const Table* old_elem_dof)
928 {
929  MFEM_VERIFY(ndofs >= old_ndofs, "Previous space is not coarser.");
930 
931  Mesh::GeometryList elem_geoms(*mesh);
932 
934  for (int i = 0; i < elem_geoms.Size(); i++)
935  {
936  GetLocalRefinementMatrices(elem_geoms[i], localP[elem_geoms[i]]);
937  }
938 
939  return RefinementMatrix_main(old_ndofs, *old_elem_dof, localP);
940 }
941 
943 (const FiniteElementSpace* fespace, Table* old_elem_dof, int old_ndofs)
944  : fespace(fespace)
945  , old_elem_dof(old_elem_dof)
946 {
947  MFEM_VERIFY(fespace->GetNDofs() >= old_ndofs,
948  "Previous space is not coarser.");
949 
950  width = old_ndofs * fespace->GetVDim();
951  height = fespace->GetVSize();
952 
953  Mesh::GeometryList elem_geoms(*fespace->GetMesh());
954 
955  for (int i = 0; i < elem_geoms.Size(); i++)
956  {
957  fespace->GetLocalRefinementMatrices(elem_geoms[i], localP[elem_geoms[i]]);
958  }
959 }
960 
962  const FiniteElementSpace *fespace, const FiniteElementSpace *coarse_fes)
963  : Operator(fespace->GetVSize(), coarse_fes->GetVSize()),
964  fespace(fespace), old_elem_dof(NULL)
965 {
966  Mesh::GeometryList elem_geoms(*fespace->GetMesh());
967 
968  for (int i = 0; i < elem_geoms.Size(); i++)
969  {
970  fespace->GetLocalRefinementMatrices(*coarse_fes, elem_geoms[i],
971  localP[elem_geoms[i]]);
972  }
973 
974  // Make a copy of the coarse elem_dof Table.
975  old_elem_dof = new Table(coarse_fes->GetElementToDofTable());
976 }
977 
979 {
980  delete old_elem_dof;
981 }
982 
984 ::Mult(const Vector &x, Vector &y) const
985 {
986  Mesh* mesh = fespace->GetMesh();
987  const CoarseFineTransformations &rtrans = mesh->GetRefinementTransforms();
988 
989  Array<int> dofs, old_dofs, old_vdofs;
990 
991  Array<char> processed(fespace->GetVSize());
992  processed = 0;
993 
994  int vdim = fespace->GetVDim();
995  int old_ndofs = width / vdim;
996 
997  for (int k = 0; k < mesh->GetNE(); k++)
998  {
999  const Embedding &emb = rtrans.embeddings[k];
1000  const Geometry::Type geom = mesh->GetElementBaseGeometry(k);
1001  const DenseMatrix &lP = localP[geom](emb.matrix);
1002 
1003  fespace->GetElementDofs(k, dofs);
1004  old_elem_dof->GetRow(emb.parent, old_dofs);
1005 
1006  for (int vd = 0; vd < vdim; vd++)
1007  {
1008  old_dofs.Copy(old_vdofs);
1009  fespace->DofsToVDofs(vd, old_vdofs, old_ndofs);
1010 
1011  for (int i = 0; i < dofs.Size(); i++)
1012  {
1013  double rsign, osign;
1014  int r = fespace->DofToVDof(dofs[i], vd);
1015  r = DecodeDof(r, rsign);
1016 
1017  if (!processed[r])
1018  {
1019  double value = 0.0;
1020  for (int j = 0; j < old_vdofs.Size(); j++)
1021  {
1022  int o = DecodeDof(old_vdofs[j], osign);
1023  value += x[o] * lP(i, j) * osign;
1024  }
1025  y[r] = value * rsign;
1026  processed[r] = 1;
1027  }
1028  }
1029  }
1030  }
1031 }
1032 
1034  const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes,
1035  BilinearFormIntegrator *mass_integ)
1036  : Operator(c_fes->GetVSize(), f_fes->GetVSize()),
1037  fine_fes(f_fes)
1038 {
1039  MFEM_VERIFY(c_fes->GetOrdering() == f_fes->GetOrdering() &&
1040  c_fes->GetVDim() == f_fes->GetVDim(),
1041  "incompatible coarse and fine FE spaces");
1042 
1044  Mesh *f_mesh = f_fes->GetMesh();
1045  const CoarseFineTransformations &rtrans = f_mesh->GetRefinementTransforms();
1046 
1047  Mesh::GeometryList elem_geoms(*f_mesh);
1049  for (int gi = 0; gi < elem_geoms.Size(); gi++)
1050  {
1051  const Geometry::Type geom = elem_geoms[gi];
1052  DenseTensor &lP = localP[geom], &lM = localM[geom];
1053  const FiniteElement *fine_fe =
1054  f_fes->fec->FiniteElementForGeometry(geom);
1055  const FiniteElement *coarse_fe =
1056  c_fes->fec->FiniteElementForGeometry(geom);
1057  const DenseTensor &pmats = rtrans.GetPointMatrices(geom);
1058 
1059  lP.SetSize(fine_fe->GetDof(), coarse_fe->GetDof(), pmats.SizeK());
1060  lM.SetSize(fine_fe->GetDof(), fine_fe->GetDof(), pmats.SizeK());
1061  emb_tr.SetIdentityTransformation(geom);
1062  for (int i = 0; i < pmats.SizeK(); i++)
1063  {
1064  emb_tr.GetPointMat() = pmats(i);
1065  emb_tr.FinalizeTransformation();
1066  // Get the local interpolation matrix for this refinement type
1067  fine_fe->GetTransferMatrix(*coarse_fe, emb_tr, lP(i));
1068  // Get the local mass matrix for this refinement type
1069  mass_integ->AssembleElementMatrix(*fine_fe, emb_tr, lM(i));
1070  }
1071  }
1072 
1073  Table ref_type_to_matrix;
1074  rtrans.GetCoarseToFineMap(*f_mesh, coarse_to_fine, coarse_to_ref_type,
1075  ref_type_to_matrix, ref_type_to_geom);
1076  MFEM_ASSERT(coarse_to_fine.Size() == c_fes->GetNE(), "");
1077 
1078  const int total_ref_types = ref_type_to_geom.Size();
1079  int num_ref_types[Geometry::NumGeom], num_fine_elems[Geometry::NumGeom];
1080  Array<int> ref_type_to_coarse_elem_offset(total_ref_types);
1081  ref_type_to_fine_elem_offset.SetSize(total_ref_types);
1082  std::fill(num_ref_types, num_ref_types+Geometry::NumGeom, 0);
1083  std::fill(num_fine_elems, num_fine_elems+Geometry::NumGeom, 0);
1084  for (int i = 0; i < total_ref_types; i++)
1085  {
1086  Geometry::Type g = ref_type_to_geom[i];
1087  ref_type_to_coarse_elem_offset[i] = num_ref_types[g];
1088  ref_type_to_fine_elem_offset[i] = num_fine_elems[g];
1089  num_ref_types[g]++;
1090  num_fine_elems[g] += ref_type_to_matrix.RowSize(i);
1091  }
1092  DenseTensor localPtMP[Geometry::NumGeom];
1093  for (int g = 0; g < Geometry::NumGeom; g++)
1094  {
1095  if (num_ref_types[g] == 0) { continue; }
1096  const int fine_dofs = localP[g].SizeI();
1097  const int coarse_dofs = localP[g].SizeJ();
1098  localPtMP[g].SetSize(coarse_dofs, coarse_dofs, num_ref_types[g]);
1099  localR[g].SetSize(coarse_dofs, fine_dofs, num_fine_elems[g]);
1100  }
1101  for (int i = 0; i < total_ref_types; i++)
1102  {
1103  Geometry::Type g = ref_type_to_geom[i];
1104  DenseMatrix &lPtMP = localPtMP[g](ref_type_to_coarse_elem_offset[i]);
1105  int lR_offset = ref_type_to_fine_elem_offset[i]; // offset in localR[g]
1106  const int *mi = ref_type_to_matrix.GetRow(i);
1107  const int nm = ref_type_to_matrix.RowSize(i);
1108  lPtMP = 0.0;
1109  for (int s = 0; s < nm; s++)
1110  {
1111  DenseMatrix &lP = localP[g](mi[s]);
1112  DenseMatrix &lM = localM[g](mi[s]);
1113  DenseMatrix &lR = localR[g](lR_offset+s);
1114  MultAtB(lP, lM, lR); // lR = lP^T lM
1115  AddMult(lR, lP, lPtMP); // lPtMP += lP^T lM lP
1116  }
1117  DenseMatrixInverse lPtMP_inv(lPtMP);
1118  for (int s = 0; s < nm; s++)
1119  {
1120  DenseMatrix &lR = localR[g](lR_offset+s);
1121  lPtMP_inv.Mult(lR); // lR <- (P^T M P)^{-1} P^T M
1122  }
1123  }
1124 
1125  // Make a copy of the coarse element-to-dof Table.
1126  coarse_elem_dof = new Table(c_fes->GetElementToDofTable());
1127 }
1128 
1130 {
1131  delete coarse_elem_dof;
1132 }
1133 
1135 ::Mult(const Vector &x, Vector &y) const
1136 {
1137  Array<int> c_vdofs, f_vdofs;
1138  Vector loc_x, loc_y;
1139  DenseMatrix loc_x_mat, loc_y_mat;
1140  const int vdim = fine_fes->GetVDim();
1141  const int coarse_ndofs = height/vdim;
1142  for (int coarse_el = 0; coarse_el < coarse_to_fine.Size(); coarse_el++)
1143  {
1144  coarse_elem_dof->GetRow(coarse_el, c_vdofs);
1145  fine_fes->DofsToVDofs(c_vdofs, coarse_ndofs);
1146  loc_y.SetSize(c_vdofs.Size());
1147  loc_y = 0.0;
1148  loc_y_mat.UseExternalData(loc_y.GetData(), c_vdofs.Size()/vdim, vdim);
1149  const int ref_type = coarse_to_ref_type[coarse_el];
1150  const Geometry::Type geom = ref_type_to_geom[ref_type];
1151  const int *fine_elems = coarse_to_fine.GetRow(coarse_el);
1152  const int num_fine_elems = coarse_to_fine.RowSize(coarse_el);
1153  const int lR_offset = ref_type_to_fine_elem_offset[ref_type];
1154  for (int s = 0; s < num_fine_elems; s++)
1155  {
1156  const DenseMatrix &lR = localR[geom](lR_offset+s);
1157  fine_fes->GetElementVDofs(fine_elems[s], f_vdofs);
1158  x.GetSubVector(f_vdofs, loc_x);
1159  loc_x_mat.UseExternalData(loc_x.GetData(), f_vdofs.Size()/vdim, vdim);
1160  AddMult(lR, loc_x_mat, loc_y_mat);
1161  }
1162  y.SetSubVector(c_vdofs, loc_y);
1163  }
1164 }
1165 
1167  DenseTensor &localR) const
1168 {
1169  const FiniteElement *fe = fec->FiniteElementForGeometry(geom);
1170 
1171  const CoarseFineTransformations &dtrans =
1173  const DenseTensor &pmats = dtrans.GetPointMatrices(geom);
1174 
1175  const int nmat = pmats.SizeK();
1176  const int ldof = fe->GetDof();
1177 
1179  isotr.SetIdentityTransformation(geom);
1180 
1181  // calculate local restriction matrices for all refinement types
1182  localR.SetSize(ldof, ldof, nmat);
1183  for (int i = 0; i < nmat; i++)
1184  {
1185  isotr.GetPointMat() = pmats(i);
1186  isotr.FinalizeTransformation();
1187 
1188  fe->GetLocalRestriction(isotr, localR(i));
1189  }
1190 }
1191 
1193  const Table* old_elem_dof)
1194 {
1195  MFEM_VERIFY(Nonconforming(), "Not implemented for conforming meshes.");
1196  MFEM_VERIFY(old_ndofs, "Missing previous (finer) space.");
1197  MFEM_VERIFY(ndofs <= old_ndofs, "Previous space is not finer.");
1198 
1199  Array<int> dofs, old_dofs, old_vdofs;
1200  Vector row;
1201 
1202  Mesh::GeometryList elem_geoms(*mesh);
1203 
1205  for (int i = 0; i < elem_geoms.Size(); i++)
1206  {
1207  GetLocalDerefinementMatrices(elem_geoms[i], localR[elem_geoms[i]]);
1208  }
1209 
1210  SparseMatrix *R;
1211  if (elem_geoms.Size() == 1)
1212  {
1213  R = new SparseMatrix(ndofs*vdim, old_ndofs*vdim,
1214  localR[elem_geoms[0]].SizeI());
1215  }
1216  else
1217  {
1218  R = new SparseMatrix(ndofs*vdim, old_ndofs*vdim);
1219  }
1220  Array<int> mark(R->Height());
1221  mark = 0;
1222 
1223  const CoarseFineTransformations &dtrans =
1225 
1226  MFEM_ASSERT(dtrans.embeddings.Size() == old_elem_dof->Size(), "");
1227 
1228  int num_marked = 0;
1229  for (int k = 0; k < dtrans.embeddings.Size(); k++)
1230  {
1231  const Embedding &emb = dtrans.embeddings[k];
1233  DenseMatrix &lR = localR[geom](emb.matrix);
1234 
1235  elem_dof->GetRow(emb.parent, dofs);
1236  old_elem_dof->GetRow(k, old_dofs);
1237 
1238  for (int vd = 0; vd < vdim; vd++)
1239  {
1240  old_dofs.Copy(old_vdofs);
1241  DofsToVDofs(vd, old_vdofs, old_ndofs);
1242 
1243  for (int i = 0; i < lR.Height(); i++)
1244  {
1245  if (lR(i, 0) == infinity()) { continue; }
1246 
1247  int r = DofToVDof(dofs[i], vd);
1248  int m = (r >= 0) ? r : (-1 - r);
1249 
1250  if (!mark[m])
1251  {
1252  lR.GetRow(i, row);
1253  R->SetRow(r, old_vdofs, row);
1254  mark[m] = 1;
1255  num_marked++;
1256  }
1257  }
1258  }
1259  }
1260 
1261  MFEM_VERIFY(num_marked == R->Height(),
1262  "internal error: not all rows of R were set.");
1263  if (elem_geoms.Size() != 1) { R->Finalize(); }
1264  return R;
1265 }
1266 
1268  const FiniteElementSpace &coarse_fes, Geometry::Type geom,
1269  DenseTensor &localP) const
1270 {
1271  // Assumptions: see the declaration of the method.
1272 
1273  const FiniteElement *fine_fe = fec->FiniteElementForGeometry(geom);
1274  const FiniteElement *coarse_fe =
1275  coarse_fes.fec->FiniteElementForGeometry(geom);
1276 
1278  const DenseTensor &pmats = rtrans.GetPointMatrices(geom);
1279 
1280  int nmat = pmats.SizeK();
1281 
1283  isotr.SetIdentityTransformation(geom);
1284 
1285  // Calculate the local interpolation matrices for all refinement types
1286  localP.SetSize(fine_fe->GetDof(), coarse_fe->GetDof(), nmat);
1287  for (int i = 0; i < nmat; i++)
1288  {
1289  isotr.GetPointMat() = pmats(i);
1290  isotr.FinalizeTransformation();
1291  fine_fe->GetTransferMatrix(*coarse_fe, isotr, localP(i));
1292  }
1293 }
1294 
1297  int vdim, int ordering)
1298 {
1299  this->mesh = mesh;
1300  this->fec = fec;
1301  this->vdim = vdim;
1302  this->ordering = (Ordering::Type) ordering;
1303 
1304  elem_dof = NULL;
1305  sequence = mesh->GetSequence();
1307 
1308  const NURBSFECollection *nurbs_fec =
1309  dynamic_cast<const NURBSFECollection *>(fec);
1310  if (nurbs_fec)
1311  {
1312  if (!mesh->NURBSext)
1313  {
1314  mfem_error("FiniteElementSpace::FiniteElementSpace :\n"
1315  " NURBS FE space requires NURBS mesh.");
1316  }
1317 
1318  if (NURBSext == NULL)
1319  {
1320  this->NURBSext = mesh->NURBSext;
1321  own_ext = 0;
1322  }
1323  else
1324  {
1325  this->NURBSext = NURBSext;
1326  own_ext = 1;
1327  }
1328  UpdateNURBS();
1329  cP = cR = NULL;
1330  cP_is_set = false;
1331  }
1332  else
1333  {
1334  this->NURBSext = NULL;
1335  own_ext = 0;
1336  Construct();
1337  }
1339 }
1340 
1342 {
1343  if (NURBSext && !own_ext)
1344  {
1345  mfem_error("FiniteElementSpace::StealNURBSext");
1346  }
1347  own_ext = 0;
1348 
1349  return NURBSext;
1350 }
1351 
1353 {
1354  nvdofs = 0;
1355  nedofs = 0;
1356  nfdofs = 0;
1357  nbdofs = 0;
1358  fdofs = NULL;
1359  bdofs = NULL;
1360 
1361  dynamic_cast<const NURBSFECollection *>(fec)->Reset();
1362 
1363  ndofs = NURBSext->GetNDof();
1366 }
1367 
1369 {
1370  // This method should be used only for non-NURBS spaces.
1371  MFEM_ASSERT(!NURBSext, "internal error");
1372 
1373  elem_dof = NULL;
1374  bdrElem_dof = NULL;
1375 
1377 
1378  if ( mesh->Dimension() > 1 )
1379  {
1381  }
1382  else
1383  {
1384  nedofs = 0;
1385  }
1386 
1387  ndofs = 0;
1388  nfdofs = 0;
1389  nbdofs = 0;
1390  bdofs = NULL;
1391  fdofs = NULL;
1392  cP = NULL;
1393  cR = NULL;
1394  cP_is_set = false;
1395  // Th is initialized/destroyed before this method is called.
1396 
1397  if (mesh->GetNFaces() > 0)
1398  {
1399  bool have_face_dofs = false;
1400  for (int g = Geometry::DimStart[2]; g < Geometry::DimStart[3]; g++)
1401  {
1402  if (mesh->HasGeometry(Geometry::Type(g)) &&
1404  {
1405  have_face_dofs = true;
1406  break;
1407  }
1408  }
1409  if (have_face_dofs)
1410  {
1411  fdofs = new int[mesh->GetNFaces()+1];
1412  fdofs[0] = 0;
1413  for (int i = 0; i < mesh->GetNFaces(); i++)
1414  {
1416  fdofs[i+1] = nfdofs;
1417  }
1418  }
1419  }
1420 
1421  if (mesh->Dimension() > 0)
1422  {
1423  bdofs = new int[mesh->GetNE()+1];
1424  bdofs[0] = 0;
1425  for (int i = 0; i < mesh->GetNE(); i++)
1426  {
1428  nbdofs += fec->DofForGeometry(geom);
1429  bdofs[i+1] = nbdofs;
1430  }
1431  }
1432 
1433  ndofs = nvdofs + nedofs + nfdofs + nbdofs;
1434 
1435  // Do not build elem_dof Table here: in parallel it has to be constructed
1436  // later.
1437 }
1438 
1440 {
1441  if (elem_dof)
1442  {
1443  elem_dof -> GetRow (i, dofs);
1444  }
1445  else
1446  {
1447  Array<int> V, E, Eo, F, Fo;
1448  int k, j, nv, ne, nf, nb, nfd, nd, dim;
1449  const int *ind;
1450 
1451  dim = mesh->Dimension();
1453  ne = (dim > 1) ? ( fec->DofForGeometry(Geometry::SEGMENT) ) : ( 0 );
1454  nb = (dim > 0) ? fec->DofForGeometry(mesh->GetElementBaseGeometry(i)) : 0;
1455  if (nv > 0)
1456  {
1457  mesh->GetElementVertices(i, V);
1458  }
1459  if (ne > 0)
1460  {
1461  mesh->GetElementEdges(i, E, Eo);
1462  }
1463  nfd = 0;
1464  if (dim == 3)
1465  {
1467  {
1468  mesh->GetElementFaces(i, F, Fo);
1469  for (k = 0; k < F.Size(); k++)
1470  {
1471  nfd += fec->DofForGeometry(mesh->GetFaceBaseGeometry(F[k]));
1472  }
1473  }
1474  }
1475  nd = V.Size() * nv + E.Size() * ne + nfd + nb;
1476  dofs.SetSize(nd);
1477  if (nv > 0)
1478  {
1479  for (k = 0; k < V.Size(); k++)
1480  {
1481  for (j = 0; j < nv; j++)
1482  {
1483  dofs[k*nv+j] = V[k]*nv+j;
1484  }
1485  }
1486  nv *= V.Size();
1487  }
1488  if (ne > 0)
1489  {
1490  // if (dim > 1)
1491  for (k = 0; k < E.Size(); k++)
1492  {
1494  for (j = 0; j < ne; j++)
1495  {
1496  if (ind[j] < 0)
1497  {
1498  dofs[nv+k*ne+j] = -1 - ( nvdofs+E[k]*ne+(-1-ind[j]) );
1499  }
1500  else
1501  {
1502  dofs[nv+k*ne+j] = nvdofs+E[k]*ne+ind[j];
1503  }
1504  }
1505  }
1506  }
1507  ne = nv + ne * E.Size();
1508  if (nfd > 0)
1509  // if (dim == 3)
1510  {
1511  for (k = 0; k < F.Size(); k++)
1512  {
1514  Fo[k]);
1516  for (j = 0; j < nf; j++)
1517  {
1518  if (ind[j] < 0)
1519  {
1520  dofs[ne+j] = -1 - ( nvdofs+nedofs+fdofs[F[k]]+(-1-ind[j]) );
1521  }
1522  else
1523  {
1524  dofs[ne+j] = nvdofs+nedofs+fdofs[F[k]]+ind[j];
1525  }
1526  }
1527  ne += nf;
1528  }
1529  }
1530  if (nb > 0)
1531  {
1532  k = nvdofs + nedofs + nfdofs + bdofs[i];
1533  for (j = 0; j < nb; j++)
1534  {
1535  dofs[ne+j] = k + j;
1536  }
1537  }
1538  }
1539 }
1540 
1542 {
1543  if (i < 0 || !mesh->GetNE()) { return NULL; }
1544  MFEM_VERIFY(i < mesh->GetNE(),
1545  "Invalid element id " << i << ", maximum allowed " << mesh->GetNE()-1);
1546 
1547  const FiniteElement *FE =
1549 
1550  if (NURBSext)
1551  {
1552  NURBSext->LoadFE(i, FE);
1553  }
1554 
1555  return FE;
1556 }
1557 
1559 {
1560  if (bdrElem_dof)
1561  {
1562  bdrElem_dof->GetRow(i, dofs);
1563  }
1564  else
1565  {
1566  Array<int> V, E, Eo;
1567  int k, j, nv, ne, nf, nd, iF, oF, dim;
1568  const int *ind;
1569 
1570  dim = mesh->Dimension();
1572  if (nv > 0)
1573  {
1574  mesh->GetBdrElementVertices(i, V);
1575  }
1576  ne = (dim > 1) ? ( fec->DofForGeometry(Geometry::SEGMENT) ) : ( 0 );
1577  if (ne > 0)
1578  {
1579  mesh->GetBdrElementEdges(i, E, Eo);
1580  }
1581  nd = V.Size() * nv + E.Size() * ne;
1582  nf = (dim == 3) ? (fec->DofForGeometry(
1583  mesh->GetBdrElementBaseGeometry(i))) : (0);
1584  if (nf > 0)
1585  {
1586  nd += nf;
1587  mesh->GetBdrElementFace(i, &iF, &oF);
1588  }
1589  dofs.SetSize(nd);
1590  if (nv > 0)
1591  {
1592  for (k = 0; k < V.Size(); k++)
1593  {
1594  for (j = 0; j < nv; j++)
1595  {
1596  dofs[k*nv+j] = V[k]*nv+j;
1597  }
1598  }
1599  nv *= V.Size();
1600  }
1601  if (ne > 0)
1602  {
1603  // if (dim > 1)
1604  for (k = 0; k < E.Size(); k++)
1605  {
1607  for (j = 0; j < ne; j++)
1608  {
1609  if (ind[j] < 0)
1610  {
1611  dofs[nv+k*ne+j] = -1 - ( nvdofs+E[k]*ne+(-1-ind[j]) );
1612  }
1613  else
1614  {
1615  dofs[nv+k*ne+j] = nvdofs+E[k]*ne+ind[j];
1616  }
1617  }
1618  }
1619  }
1620  if (nf > 0)
1621  // if (dim == 3)
1622  {
1623  ne = nv + ne * E.Size();
1624  ind = fec->DofOrderForOrientation(
1626  for (j = 0; j < nf; j++)
1627  {
1628  if (ind[j] < 0)
1629  {
1630  dofs[ne+j] = -1 - ( nvdofs+nedofs+fdofs[iF]+(-1-ind[j]) );
1631  }
1632  else
1633  {
1634  dofs[ne+j] = nvdofs+nedofs+fdofs[iF]+ind[j];
1635  }
1636  }
1637  }
1638  }
1639 }
1640 
1642 {
1643  int j, k, nv, ne, nf, nd, dim = mesh->Dimension();
1644  Array<int> V, E, Eo;
1645  const int *ind;
1646 
1647  // for 1D, 2D and 3D faces
1649  ne = (dim > 1) ? fec->DofForGeometry(Geometry::SEGMENT) : 0;
1650  if (nv > 0)
1651  {
1652  mesh->GetFaceVertices(i, V);
1653  }
1654  if (ne > 0)
1655  {
1656  mesh->GetFaceEdges(i, E, Eo);
1657  }
1658  nf = (fdofs) ? (fdofs[i+1]-fdofs[i]) : (0);
1659  nd = V.Size() * nv + E.Size() * ne + nf;
1660  dofs.SetSize(nd);
1661  if (nv > 0)
1662  {
1663  for (k = 0; k < V.Size(); k++)
1664  {
1665  for (j = 0; j < nv; j++)
1666  {
1667  dofs[k*nv+j] = V[k]*nv+j;
1668  }
1669  }
1670  }
1671  nv *= V.Size();
1672  if (ne > 0)
1673  {
1674  for (k = 0; k < E.Size(); k++)
1675  {
1677  for (j = 0; j < ne; j++)
1678  {
1679  if (ind[j] < 0)
1680  {
1681  dofs[nv+k*ne+j] = -1 - ( nvdofs+E[k]*ne+(-1-ind[j]) );
1682  }
1683  else
1684  {
1685  dofs[nv+k*ne+j] = nvdofs+E[k]*ne+ind[j];
1686  }
1687  }
1688  }
1689  }
1690  ne = nv + ne * E.Size();
1691  if (nf > 0)
1692  {
1693  for (j = nvdofs+nedofs+fdofs[i], k = 0; k < nf; j++, k++)
1694  {
1695  dofs[ne+k] = j;
1696  }
1697  }
1698 }
1699 
1701 {
1702  int j, k, nv, ne;
1703  Array<int> V;
1704 
1706  if (nv > 0)
1707  {
1708  mesh->GetEdgeVertices(i, V);
1709  }
1711  dofs.SetSize(2*nv+ne);
1712  if (nv > 0)
1713  {
1714  for (k = 0; k < 2; k++)
1715  {
1716  for (j = 0; j < nv; j++)
1717  {
1718  dofs[k*nv+j] = V[k]*nv+j;
1719  }
1720  }
1721  }
1722  nv *= 2;
1723  for (j = 0, k = nvdofs+i*ne; j < ne; j++, k++)
1724  {
1725  dofs[nv+j] = k;
1726  }
1727 }
1728 
1730 {
1731  int j, nv;
1732 
1734  dofs.SetSize(nv);
1735  for (j = 0; j < nv; j++)
1736  {
1737  dofs[j] = i*nv+j;
1738  }
1739 }
1740 
1742 {
1743  int j, k, nb;
1744  if (mesh->Dimension() == 0) { dofs.SetSize(0); return; }
1745  nb = fec -> DofForGeometry (mesh -> GetElementBaseGeometry (i));
1746  dofs.SetSize (nb);
1747  k = nvdofs + nedofs + nfdofs + bdofs[i];
1748  for (j = 0; j < nb; j++)
1749  {
1750  dofs[j] = k + j;
1751  }
1752 }
1753 
1755 {
1756  int j, k, ne;
1757 
1758  ne = fec -> DofForGeometry (Geometry::SEGMENT);
1759  dofs.SetSize (ne);
1760  for (j = 0, k = nvdofs+i*ne; j < ne; j++, k++)
1761  {
1762  dofs[j] = k;
1763  }
1764 }
1765 
1767 {
1768  int j, k, nf;
1769 
1770  nf = (fdofs) ? (fdofs[i+1]-fdofs[i]) : (0);
1771  dofs.SetSize (nf);
1772  if (nf > 0)
1773  {
1774  for (j = 0, k = nvdofs+nedofs+fdofs[i]; j < nf; j++, k++)
1775  {
1776  dofs[j] = k;
1777  }
1778  }
1779 }
1780 
1782 {
1783  const FiniteElement *BE;
1784 
1785  switch ( mesh->Dimension() )
1786  {
1787  case 1:
1789  break;
1790  case 2:
1792  break;
1793  case 3:
1794  default:
1797  }
1798 
1799  if (NURBSext)
1800  {
1801  NURBSext->LoadBE(i, BE);
1802  }
1803 
1804  return BE;
1805 }
1806 
1808 {
1809  const FiniteElement *fe;
1810 
1811  switch (mesh->Dimension())
1812  {
1813  case 1:
1815  break;
1816  case 2:
1818  break;
1819  case 3:
1820  default:
1822  }
1823 
1824  // if (NURBSext)
1825  // NURBSext->LoadFaceElement(i, fe);
1826 
1827  return fe;
1828 }
1829 
1831 {
1833 }
1834 
1836  int i, Geometry::Type geom_type) const
1837 {
1838  return fec->TraceFiniteElementForGeometry(geom_type);
1839 }
1840 
1842 {
1843  Destroy();
1844 }
1845 
1847 {
1848  delete cR;
1849  delete cP;
1850  Th.Clear();
1851  L2E_nat.Clear();
1852  L2E_lex.Clear();
1853  for (int i = 0; i < E2Q_array.Size(); i++)
1854  {
1855  delete E2Q_array[i];
1856  }
1857  E2Q_array.SetSize(0);
1858 
1861 
1862  if (NURBSext)
1863  {
1864  if (own_ext) { delete NURBSext; }
1865  }
1866  else
1867  {
1868  delete elem_dof;
1869  delete bdrElem_dof;
1870 
1871  delete [] bdofs;
1872  delete [] fdofs;
1873  }
1874 }
1875 
1877  const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
1878 {
1879  // Assumptions: see the declaration of the method.
1880 
1881  if (T.Type() == Operator::MFEM_SPARSEMAT)
1882  {
1883  Mesh::GeometryList elem_geoms(*mesh);
1884 
1886  for (int i = 0; i < elem_geoms.Size(); i++)
1887  {
1888  GetLocalRefinementMatrices(coarse_fes, elem_geoms[i],
1889  localP[elem_geoms[i]]);
1890  }
1891  T.Reset(RefinementMatrix_main(coarse_fes.GetNDofs(),
1892  coarse_fes.GetElementToDofTable(),
1893  localP));
1894  }
1895  else
1896  {
1897  T.Reset(new RefinementOperator(this, &coarse_fes));
1898  }
1899 }
1900 
1902  const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
1903 {
1904  const SparseMatrix *coarse_P = coarse_fes.GetConformingProlongation();
1905 
1906  Operator::Type req_type = T.Type();
1907  GetTransferOperator(coarse_fes, T);
1908 
1909  if (req_type == Operator::MFEM_SPARSEMAT)
1910  {
1912  {
1913  T.Reset(mfem::Mult(*cR, *T.As<SparseMatrix>()));
1914  }
1915  if (coarse_P)
1916  {
1917  T.Reset(mfem::Mult(*T.As<SparseMatrix>(), *coarse_P));
1918  }
1919  }
1920  else
1921  {
1922  const int RP_case = bool(GetConformingRestriction()) + 2*bool(coarse_P);
1923  if (RP_case == 0) { return; }
1924  const bool owner = T.OwnsOperator();
1925  T.SetOperatorOwner(false);
1926  switch (RP_case)
1927  {
1928  case 1:
1929  T.Reset(new ProductOperator(cR, T.Ptr(), false, owner));
1930  break;
1931  case 2:
1932  T.Reset(new ProductOperator(T.Ptr(), coarse_P, owner, false));
1933  break;
1934  case 3:
1936  cR, T.Ptr(), coarse_P, false, owner, false));
1937  break;
1938  }
1939  }
1940 }
1941 
1942 void FiniteElementSpace::Update(bool want_transform)
1943 {
1944  if (mesh->GetSequence() == sequence)
1945  {
1946  return; // mesh and space are in sync, no-op
1947  }
1948  if (want_transform && mesh->GetSequence() != sequence + 1)
1949  {
1950  MFEM_ABORT("Error in update sequence. Space needs to be updated after "
1951  "each mesh modification.");
1952  }
1953  sequence = mesh->GetSequence();
1954 
1955  if (NURBSext)
1956  {
1957  UpdateNURBS();
1958  return;
1959  }
1960 
1961  Table* old_elem_dof = NULL;
1962  int old_ndofs;
1963 
1964  // save old DOF table
1965  if (want_transform)
1966  {
1967  old_elem_dof = elem_dof;
1968  elem_dof = NULL;
1969  old_ndofs = ndofs;
1970  }
1971 
1972  Destroy(); // calls Th.Clear()
1973  Construct();
1975 
1976  if (want_transform)
1977  {
1978  // calculate appropriate GridFunction transformation
1979  switch (mesh->GetLastOperation())
1980  {
1981  case Mesh::REFINE:
1982  {
1983  if (Th.Type() != Operator::MFEM_SPARSEMAT)
1984  {
1985  Th.Reset(new RefinementOperator(this, old_elem_dof, old_ndofs));
1986  // The RefinementOperator takes ownership of 'old_elem_dof', so
1987  // we no longer own it:
1988  old_elem_dof = NULL;
1989  }
1990  else
1991  {
1992  // calculate fully assembled matrix
1993  Th.Reset(RefinementMatrix(old_ndofs, old_elem_dof));
1994  }
1995  break;
1996  }
1997 
1998  case Mesh::DEREFINE:
1999  {
2001  Th.Reset(DerefinementMatrix(old_ndofs, old_elem_dof));
2002  if (cP && cR)
2003  {
2004  Th.SetOperatorOwner(false);
2006  false, false, true));
2007  }
2008  break;
2009  }
2010 
2011  default:
2012  break;
2013  }
2014 
2015  delete old_elem_dof;
2016  }
2017 }
2018 
2019 void FiniteElementSpace::Save(std::ostream &out) const
2020 {
2021  int fes_format = 90; // the original format, v0.9
2022  bool nurbs_unit_weights = false;
2023 
2024  // Determine the format that should be used.
2025  if (!NURBSext)
2026  {
2027  // TODO: if this is a variable-order FE space, use fes_format = 100.
2028  }
2029  else
2030  {
2031  const NURBSFECollection *nurbs_fec =
2032  dynamic_cast<const NURBSFECollection *>(fec);
2033  MFEM_VERIFY(nurbs_fec, "invalid FE collection");
2034  nurbs_fec->SetOrder(NURBSext->GetOrder());
2035  const double eps = 5e-14;
2036  nurbs_unit_weights = (NURBSext->GetWeights().Min() >= 1.0-eps &&
2037  NURBSext->GetWeights().Max() <= 1.0+eps);
2039  (NURBSext != mesh->NURBSext && !nurbs_unit_weights) ||
2040  (NURBSext->GetMaster().Size() != 0 ))
2041  {
2042  fes_format = 100; // v1.0 format
2043  }
2044  }
2045 
2046  out << (fes_format == 90 ?
2047  "FiniteElementSpace\n" : "MFEM FiniteElementSpace v1.0\n")
2048  << "FiniteElementCollection: " << fec->Name() << '\n'
2049  << "VDim: " << vdim << '\n'
2050  << "Ordering: " << ordering << '\n';
2051 
2052  if (fes_format == 100) // v1.0
2053  {
2054  if (!NURBSext)
2055  {
2056  // TODO: this is a variable-order FE space --> write 'element_orders'.
2057  }
2058  else if (NURBSext != mesh->NURBSext)
2059  {
2061  {
2062  out << "NURBS_order\n" << NURBSext->GetOrder() << '\n';
2063  }
2064  else
2065  {
2066  out << "NURBS_orders\n";
2067  // 1 = do not write the size, just the entries:
2068  NURBSext->GetOrders().Save(out, 1);
2069  }
2070  // If periodic BCs are given, write connectivity
2071  if (NURBSext->GetMaster().Size() != 0 )
2072  {
2073  out <<"NURBS_periodic\n";
2074  NURBSext->GetMaster().Save(out);
2075  NURBSext->GetSlave().Save(out);
2076  }
2077  // If the weights are not unit, write them to the output:
2078  if (!nurbs_unit_weights)
2079  {
2080  out << "NURBS_weights\n";
2081  NURBSext->GetWeights().Print(out, 1);
2082  }
2083  }
2084  out << "End: MFEM FiniteElementSpace v1.0\n";
2085  }
2086 }
2087 
2089 {
2090  string buff;
2091  int fes_format = 0, ord;
2092  FiniteElementCollection *r_fec;
2093 
2094  Destroy();
2095 
2096  input >> std::ws;
2097  getline(input, buff); // 'FiniteElementSpace'
2098  filter_dos(buff);
2099  if (buff == "FiniteElementSpace") { fes_format = 90; /* v0.9 */ }
2100  else if (buff == "MFEM FiniteElementSpace v1.0") { fes_format = 100; }
2101  else { MFEM_ABORT("input stream is not a FiniteElementSpace!"); }
2102  getline(input, buff, ' '); // 'FiniteElementCollection:'
2103  input >> std::ws;
2104  getline(input, buff);
2105  filter_dos(buff);
2106  r_fec = FiniteElementCollection::New(buff.c_str());
2107  getline(input, buff, ' '); // 'VDim:'
2108  input >> vdim;
2109  getline(input, buff, ' '); // 'Ordering:'
2110  input >> ord;
2111 
2112  NURBSFECollection *nurbs_fec = dynamic_cast<NURBSFECollection*>(r_fec);
2113  NURBSExtension *NURBSext = NULL;
2114  if (fes_format == 90) // original format, v0.9
2115  {
2116  if (nurbs_fec)
2117  {
2118  MFEM_VERIFY(m->NURBSext, "NURBS FE collection requires a NURBS mesh!");
2119  const int order = nurbs_fec->GetOrder();
2120  if (order != m->NURBSext->GetOrder() &&
2122  {
2123  NURBSext = new NURBSExtension(m->NURBSext, order);
2124  }
2125  }
2126  }
2127  else if (fes_format == 100) // v1.0
2128  {
2129  while (1)
2130  {
2131  skip_comment_lines(input, '#');
2132  MFEM_VERIFY(input.good(), "error reading FiniteElementSpace v1.0");
2133  getline(input, buff);
2134  filter_dos(buff);
2135  if (buff == "NURBS_order" || buff == "NURBS_orders")
2136  {
2137  MFEM_VERIFY(nurbs_fec,
2138  buff << ": NURBS FE collection is required!");
2139  MFEM_VERIFY(m->NURBSext, buff << ": NURBS mesh is required!");
2140  MFEM_VERIFY(!NURBSext, buff << ": order redefinition!");
2141  if (buff == "NURBS_order")
2142  {
2143  int order;
2144  input >> order;
2145  NURBSext = new NURBSExtension(m->NURBSext, order);
2146  }
2147  else
2148  {
2149  Array<int> orders;
2150  orders.Load(m->NURBSext->GetNKV(), input);
2151  NURBSext = new NURBSExtension(m->NURBSext, orders);
2152  }
2153  }
2154  else if (buff == "NURBS_periodic")
2155  {
2156  Array<int> master, slave;
2157  master.Load(input);
2158  slave.Load(input);
2159  NURBSext->ConnectBoundaries(master,slave);
2160  }
2161  else if (buff == "NURBS_weights")
2162  {
2163  MFEM_VERIFY(NURBSext, "NURBS_weights: NURBS_orders have to be "
2164  "specified before NURBS_weights!");
2165  NURBSext->GetWeights().Load(input, NURBSext->GetNDof());
2166  }
2167  else if (buff == "element_orders")
2168  {
2169  MFEM_VERIFY(!nurbs_fec, "section element_orders cannot be used "
2170  "with a NURBS FE collection");
2171  MFEM_ABORT("element_orders: not implemented yet!");
2172  }
2173  else if (buff == "End: MFEM FiniteElementSpace v1.0")
2174  {
2175  break;
2176  }
2177  else
2178  {
2179  MFEM_ABORT("unknown section: " << buff);
2180  }
2181  }
2182  }
2183 
2184  Constructor(m, NURBSext, r_fec, vdim, ord);
2185 
2186  return r_fec;
2187 }
2188 
2189 
2191 {
2192  // protected method
2193  int offset = 0;
2194  const int num_elem = mesh->GetNE();
2195  element_offsets = new int[num_elem + 1];
2196  for (int g = 0; g < Geometry::NumGeom; g++)
2197  {
2198  int_rule[g] = NULL;
2199  }
2200  for (int i = 0; i < num_elem; i++)
2201  {
2202  element_offsets[i] = offset;
2203  int geom = mesh->GetElementBaseGeometry(i);
2204  if (int_rule[geom] == NULL)
2205  {
2206  int_rule[geom] = &IntRules.Get(geom, order);
2207  }
2208  offset += int_rule[geom]->GetNPoints();
2209  }
2210  element_offsets[num_elem] = size = offset;
2211 }
2212 
2213 QuadratureSpace::QuadratureSpace(Mesh *mesh_, std::istream &in)
2214  : mesh(mesh_)
2215 {
2216  const char *msg = "invalid input stream";
2217  string ident;
2218 
2219  in >> ident; MFEM_VERIFY(ident == "QuadratureSpace", msg);
2220  in >> ident; MFEM_VERIFY(ident == "Type:", msg);
2221  in >> ident;
2222  if (ident == "default_quadrature")
2223  {
2224  in >> ident; MFEM_VERIFY(ident == "Order:", msg);
2225  in >> order;
2226  }
2227  else
2228  {
2229  MFEM_ABORT("unknown QuadratureSpace type: " << ident);
2230  return;
2231  }
2232 
2233  Construct();
2234 }
2235 
2236 void QuadratureSpace::Save(std::ostream &out) const
2237 {
2238  out << "QuadratureSpace\n"
2239  << "Type: default_quadrature\n"
2240  << "Order: " << order << '\n';
2241 }
2242 
2243 
2245  FiniteElementSpace &ran_fes_)
2246  : dom_fes(dom_fes_), ran_fes(ran_fes_),
2247  oper_type(Operator::ANY_TYPE),
2248  fw_t_oper(), bw_t_oper()
2249 {
2250 #ifdef MFEM_USE_MPI
2251  const bool par_dom = dynamic_cast<ParFiniteElementSpace*>(&dom_fes);
2252  const bool par_ran = dynamic_cast<ParFiniteElementSpace*>(&ran_fes);
2253  MFEM_VERIFY(par_dom == par_ran, "the domain and range FE spaces must both"
2254  " be either serial or parallel");
2255  parallel = par_dom;
2256 #endif
2257 }
2258 
2260  FiniteElementSpace &fes_in, FiniteElementSpace &fes_out,
2261  const Operator &oper, OperatorHandle &t_oper)
2262 {
2263  if (t_oper.Ptr())
2264  {
2265  return *t_oper.Ptr();
2266  }
2267 
2268  if (!Parallel())
2269  {
2270  const SparseMatrix *in_cP = fes_in.GetConformingProlongation();
2271  const SparseMatrix *out_cR = fes_out.GetConformingRestriction();
2273  {
2274  const SparseMatrix *mat = dynamic_cast<const SparseMatrix *>(&oper);
2275  MFEM_VERIFY(mat != NULL, "Operator is not a SparseMatrix");
2276  if (!out_cR)
2277  {
2278  t_oper.Reset(const_cast<SparseMatrix*>(mat), false);
2279  }
2280  else
2281  {
2282  t_oper.Reset(mfem::Mult(*out_cR, *mat));
2283  }
2284  if (in_cP)
2285  {
2286  t_oper.Reset(mfem::Mult(*t_oper.As<SparseMatrix>(), *in_cP));
2287  }
2288  }
2289  else if (oper_type == Operator::ANY_TYPE)
2290  {
2291  const int RP_case = bool(out_cR) + 2*bool(in_cP);
2292  switch (RP_case)
2293  {
2294  case 0:
2295  t_oper.Reset(const_cast<Operator*>(&oper), false);
2296  break;
2297  case 1:
2298  t_oper.Reset(
2299  new ProductOperator(out_cR, &oper, false, false));
2300  break;
2301  case 2:
2302  t_oper.Reset(
2303  new ProductOperator(&oper, in_cP, false, false));
2304  break;
2305  case 3:
2306  t_oper.Reset(
2308  out_cR, &oper, in_cP, false, false, false));
2309  break;
2310  }
2311  }
2312  else
2313  {
2314  MFEM_ABORT("Operator::Type is not supported: " << oper_type);
2315  }
2316  }
2317  else // Parallel() == true
2318  {
2319 #ifdef MFEM_USE_MPI
2320  const SparseMatrix *out_R = fes_out.GetRestrictionMatrix();
2322  {
2323  const ParFiniteElementSpace *pfes_in =
2324  dynamic_cast<const ParFiniteElementSpace *>(&fes_in);
2325  const ParFiniteElementSpace *pfes_out =
2326  dynamic_cast<const ParFiniteElementSpace *>(&fes_out);
2327  const SparseMatrix *sp_mat = dynamic_cast<const SparseMatrix *>(&oper);
2328  const HypreParMatrix *hy_mat;
2329  if (sp_mat)
2330  {
2331  SparseMatrix *RA = mfem::Mult(*out_R, *sp_mat);
2332  t_oper.Reset(pfes_in->Dof_TrueDof_Matrix()->
2333  LeftDiagMult(*RA, pfes_out->GetTrueDofOffsets()));
2334  delete RA;
2335  }
2336  else if ((hy_mat = dynamic_cast<const HypreParMatrix *>(&oper)))
2337  {
2338  HypreParMatrix *RA =
2339  hy_mat->LeftDiagMult(*out_R, pfes_out->GetTrueDofOffsets());
2340  t_oper.Reset(mfem::ParMult(RA, pfes_in->Dof_TrueDof_Matrix()));
2341  delete RA;
2342  }
2343  else
2344  {
2345  MFEM_ABORT("unknown Operator type");
2346  }
2347  }
2348  else if (oper_type == Operator::ANY_TYPE)
2349  {
2350  t_oper.Reset(new TripleProductOperator(
2351  out_R, &oper, fes_in.GetProlongationMatrix(),
2352  false, false, false));
2353  }
2354  else
2355  {
2356  MFEM_ABORT("Operator::Type is not supported: " << oper_type);
2357  }
2358 #endif
2359  }
2360 
2361  return *t_oper.Ptr();
2362 }
2363 
2364 
2366 {
2367  if (own_mass_integ) { delete mass_integ; }
2368 }
2369 
2371  BilinearFormIntegrator *mass_integ_, bool own_mass_integ_)
2372 {
2373  if (own_mass_integ) { delete mass_integ; }
2374 
2375  mass_integ = mass_integ_;
2376  own_mass_integ = own_mass_integ_;
2377 }
2378 
2380 {
2381  if (F.Ptr())
2382  {
2383  return *F.Ptr();
2384  }
2385 
2386  // Costruct F
2388  {
2390  }
2391  else if (oper_type == Operator::MFEM_SPARSEMAT)
2392  {
2393  Mesh::GeometryList elem_geoms(*ran_fes.GetMesh());
2394 
2396  for (int i = 0; i < elem_geoms.Size(); i++)
2397  {
2399  localP[elem_geoms[i]]);
2400  }
2403  }
2404  else
2405  {
2406  MFEM_ABORT("Operator::Type is not supported: " << oper_type);
2407  }
2408 
2409  return *F.Ptr();
2410 }
2411 
2413 {
2414  if (B.Ptr())
2415  {
2416  return *B.Ptr();
2417  }
2418 
2419  // Construct B, if not set, define a suitable mass_integ
2420  if (!mass_integ && ran_fes.GetNE() > 0)
2421  {
2422  const FiniteElement *f_fe_0 = ran_fes.GetFE(0);
2423  const int map_type = f_fe_0->GetMapType();
2424  if (map_type == FiniteElement::VALUE ||
2425  map_type == FiniteElement::INTEGRAL)
2426  {
2427  mass_integ = new MassIntegrator;
2428  }
2429  else if (map_type == FiniteElement::H_DIV ||
2430  map_type == FiniteElement::H_CURL)
2431  {
2433  }
2434  else
2435  {
2436  MFEM_ABORT("unknown type of FE space");
2437  }
2438  own_mass_integ = true;
2439  }
2441  {
2443  &ran_fes, &dom_fes, mass_integ));
2444  }
2445  else
2446  {
2447  MFEM_ABORT("Operator::Type is not supported: " << oper_type);
2448  }
2449 
2450  return *B.Ptr();
2451 }
2452 
2453 
2455  const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
2456  : fes_ho(fes_ho_), fes_lor(fes_lor_)
2457 {
2458  Mesh *mesh_ho = fes_ho.GetMesh();
2459  MFEM_VERIFY(mesh_ho->GetNumGeometries(mesh_ho->Dimension()) <= 1,
2460  "mixed meshes are not supported");
2461 
2462  // If the local mesh is empty, skip all computations
2463  if (mesh_ho->GetNE() == 0) { return; }
2464 
2465  const FiniteElement *fe_lor = fes_lor.GetFE(0);
2466  const FiniteElement *fe_ho = fes_ho.GetFE(0);
2467  ndof_lor = fe_lor->GetDof();
2468  ndof_ho = fe_ho->GetDof();
2469 
2470  const int nel_lor = fes_lor.GetNE();
2471  const int nel_ho = fes_ho.GetNE();
2472 
2473  nref = nel_lor/nel_ho;
2474 
2475  // Construct the mapping from HO to LOR
2476  // ho2lor.GetRow(iho) will give all the LOR elements contained in iho
2477  ho2lor.SetSize(nel_ho, nref);
2478  const CoarseFineTransformations &cf_tr =
2479  fes_lor.GetMesh()->GetRefinementTransforms();
2480  for (int ilor=0; ilor<nel_lor; ++ilor)
2481  {
2482  int iho = cf_tr.embeddings[ilor].parent;
2483  ho2lor.AddConnection(iho, ilor);
2484  }
2485  ho2lor.ShiftUpI();
2486 
2487  // R will contain the restriction (L^2 projection operator) defined on
2488  // each coarse HO element (and corresponding patch of LOR elements)
2489  R.SetSize(ndof_lor*nref, ndof_ho, nel_ho);
2490  // P will contain the corresponding prolongation operator
2491  P.SetSize(ndof_ho, ndof_lor*nref, nel_ho);
2492 
2493  DenseMatrix Minv_lor(ndof_lor*nref, ndof_lor*nref);
2494  DenseMatrix M_mixed(ndof_lor*nref, ndof_ho);
2495 
2496  MassIntegrator mi;
2497  DenseMatrix M_lor_el(ndof_lor, ndof_lor);
2498  DenseMatrixInverse Minv_lor_el(&M_lor_el);
2499  DenseMatrix M_lor(ndof_lor*nref, ndof_lor*nref);
2500  DenseMatrix M_mixed_el(ndof_lor, ndof_ho);
2501 
2502  Minv_lor = 0.0;
2503  M_lor = 0.0;
2504 
2505  DenseMatrix RtMlor(ndof_ho, ndof_lor*nref);
2506  DenseMatrix RtMlorR(ndof_ho, ndof_ho);
2507  DenseMatrixInverse RtMlorR_inv(&RtMlorR);
2508 
2510  IsoparametricTransformation &emb_tr = ip_tr.Transf;
2511 
2512  Vector shape_ho(ndof_ho);
2513  Vector shape_lor(ndof_lor);
2514 
2515  const Geometry::Type geom = fe_ho->GetGeomType();
2516  const DenseTensor &pmats = cf_tr.GetPointMatrices(geom);
2517  emb_tr.SetIdentityTransformation(geom);
2518 
2519  for (int iho=0; iho<nel_ho; ++iho)
2520  {
2521  for (int iref=0; iref<nref; ++iref)
2522  {
2523  // Assemble the low-order refined mass matrix and invert locally
2524  int ilor = ho2lor.GetRow(iho)[iref];
2525  ElementTransformation *el_tr = fes_lor.GetElementTransformation(ilor);
2526  mi.AssembleElementMatrix(*fe_lor, *el_tr, M_lor_el);
2527  M_lor.CopyMN(M_lor_el, iref*ndof_lor, iref*ndof_lor);
2528  Minv_lor_el.Factor();
2529  Minv_lor_el.GetInverseMatrix(M_lor_el);
2530  // Insert into the diagonal of the patch LOR mass matrix
2531  Minv_lor.CopyMN(M_lor_el, iref*ndof_lor, iref*ndof_lor);
2532 
2533  // Now assemble the block-row of the mixed mass matrix associated
2534  // with integrating HO functions against LOR functions on the LOR
2535  // sub-element.
2536 
2537  // Create the transformation that embeds the fine low-order element
2538  // within the coarse high-order element in reference space
2539  emb_tr.GetPointMat() = pmats(iref);
2540  emb_tr.FinalizeTransformation();
2541 
2542  int order = fe_lor->GetOrder() + fe_ho->GetOrder() + el_tr->OrderW();
2543  const IntegrationRule *ir = &IntRules.Get(geom, order);
2544  M_mixed_el = 0.0;
2545  for (int i = 0; i < ir->GetNPoints(); i++)
2546  {
2547  const IntegrationPoint &ip_lor = ir->IntPoint(i);
2548  IntegrationPoint ip_ho;
2549  ip_tr.Transform(ip_lor, ip_ho);
2550  fe_lor->CalcShape(ip_lor, shape_lor);
2551  fe_ho->CalcShape(ip_ho, shape_ho);
2552  el_tr->SetIntPoint(&ip_lor);
2553  // For now we use the geometry information from the LOR space
2554  // which means we won't be mass conservative if the mesh is curved
2555  double w = el_tr->Weight()*ip_lor.weight;
2556  shape_lor *= w;
2557  AddMultVWt(shape_lor, shape_ho, M_mixed_el);
2558  }
2559  M_mixed.CopyMN(M_mixed_el, iref*ndof_lor, 0);
2560  }
2561  mfem::Mult(Minv_lor, M_mixed, R(iho));
2562 
2563  mfem::MultAtB(R(iho), M_lor, RtMlor);
2564  mfem::Mult(RtMlor, R(iho), RtMlorR);
2565  RtMlorR_inv.Factor();
2566  RtMlorR_inv.Mult(RtMlor, P(iho));
2567  }
2568 }
2569 
2571  const Vector &x, Vector &y) const
2572 {
2573  int vdim = fes_ho.GetVDim();
2574  Array<int> vdofs;
2575  DenseMatrix xel_mat(ndof_ho, vdim);
2576  DenseMatrix yel_mat(ndof_lor*nref, vdim);
2577  for (int iho=0; iho<fes_ho.GetNE(); ++iho)
2578  {
2579  fes_ho.GetElementVDofs(iho, vdofs);
2580  x.GetSubVector(vdofs, xel_mat.GetData());
2581  mfem::Mult(R(iho), xel_mat, yel_mat);
2582  // Place result correctly into the low-order vector
2583  for (int iref=0; iref<nref; ++iref)
2584  {
2585  int ilor = ho2lor.GetRow(iho)[iref];
2586  for (int vd=0; vd<vdim; ++vd)
2587  {
2588  fes_lor.GetElementDofs(ilor, vdofs);
2589  fes_lor.DofsToVDofs(vd, vdofs);
2590  y.SetSubVector(vdofs, &yel_mat(iref*ndof_lor,vd));
2591  }
2592  }
2593  }
2594 }
2595 
2597  const Vector &x, Vector &y) const
2598 {
2599  int vdim = fes_ho.GetVDim();
2600  Array<int> vdofs;
2601  DenseMatrix xel_mat(ndof_lor*nref, vdim);
2602  DenseMatrix yel_mat(ndof_ho, vdim);
2603  for (int iho=0; iho<fes_ho.GetNE(); ++iho)
2604  {
2605  // Extract the LOR DOFs
2606  for (int iref=0; iref<nref; ++iref)
2607  {
2608  int ilor = ho2lor.GetRow(iho)[iref];
2609  for (int vd=0; vd<vdim; ++vd)
2610  {
2611  fes_lor.GetElementDofs(ilor, vdofs);
2612  fes_lor.DofsToVDofs(vd, vdofs);
2613  x.GetSubVector(vdofs, &xel_mat(iref*ndof_lor, vd));
2614  }
2615  }
2616  // Locally prolongate
2617  mfem::Mult(P(iho), xel_mat, yel_mat);
2618  // Place the result in the HO vector
2619  fes_ho.GetElementVDofs(iho, vdofs);
2620  y.SetSubVector(vdofs, yel_mat.GetData());
2621  }
2622 }
2623 
2625 {
2626  if (!F) { F = new L2Projection(dom_fes, ran_fes); }
2627  return *F;
2628 }
2629 
2631 {
2632  if (!B)
2633  {
2634  if (!F) { F = new L2Projection(dom_fes, ran_fes); }
2635  B = new L2Prolongation(*F);
2636  }
2637  return *B;
2638 }
2639 
2640 
2642  ElementDofOrdering e_ordering)
2643  : fes(f),
2644  ne(fes.GetNE()),
2645  vdim(fes.GetVDim()),
2646  byvdim(fes.GetOrdering() == Ordering::byVDIM),
2647  ndofs(fes.GetNDofs()),
2648  dof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
2649  nedofs(ne*dof),
2650  offsets(ndofs+1),
2651  indices(ne*dof)
2652 {
2653  // Assuming all finite elements are the same.
2654  height = vdim*ne*dof;
2655  width = fes.GetVSize();
2656  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
2657  const int *dof_map = NULL;
2658  if (dof_reorder && ne > 0)
2659  {
2660  for (int e = 0; e < ne; ++e)
2661  {
2662  const FiniteElement *fe = fes.GetFE(e);
2663  const TensorBasisElement* el =
2664  dynamic_cast<const TensorBasisElement*>(fe);
2665  if (el) { continue; }
2666  mfem_error("Finite element not suitable for lexicographic ordering");
2667  }
2668  const FiniteElement *fe = fes.GetFE(0);
2669  const TensorBasisElement* el =
2670  dynamic_cast<const TensorBasisElement*>(fe);
2671  const Array<int> &fe_dof_map = el->GetDofMap();
2672  MFEM_VERIFY(fe_dof_map.Size() > 0, "invalid dof map");
2673  dof_map = fe_dof_map.GetData();
2674  }
2675  const Table& e2dTable = fes.GetElementToDofTable();
2676  const int* elementMap = e2dTable.GetJ();
2677  // We will be keeping a count of how many local nodes point to its global dof
2678  for (int i = 0; i <= ndofs; ++i)
2679  {
2680  offsets[i] = 0;
2681  }
2682  for (int e = 0; e < ne; ++e)
2683  {
2684  for (int d = 0; d < dof; ++d)
2685  {
2686  const int gid = elementMap[dof*e + d];
2687  ++offsets[gid + 1];
2688  }
2689  }
2690  // Aggregate to find offsets for each global dof
2691  for (int i = 1; i <= ndofs; ++i)
2692  {
2693  offsets[i] += offsets[i - 1];
2694  }
2695  // For each global dof, fill in all local nodes that point to it
2696  for (int e = 0; e < ne; ++e)
2697  {
2698  for (int d = 0; d < dof; ++d)
2699  {
2700  const int did = (!dof_reorder)?d:dof_map[d];
2701  const int gid = elementMap[dof*e + did];
2702  const int lid = dof*e + d;
2703  indices[offsets[gid]++] = lid;
2704  }
2705  }
2706  // We shifted the offsets vector by 1 by using it as a counter.
2707  // Now we shift it back.
2708  for (int i = ndofs; i > 0; --i)
2709  {
2710  offsets[i] = offsets[i - 1];
2711  }
2712  offsets[0] = 0;
2713 }
2714 
2715 void ElementRestriction::Mult(const Vector& x, Vector& y) const
2716 {
2717  // Assumes all elements have the same number of dofs
2718  const int nd = dof;
2719  const int vd = vdim;
2720  const bool t = byvdim;
2721  auto d_offsets = offsets.Read();
2722  auto d_indices = indices.Read();
2723  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
2724  auto d_y = Reshape(y.Write(), nd, vd, ne);
2725  MFEM_FORALL(i, ndofs,
2726  {
2727  const int offset = d_offsets[i];
2728  const int nextOffset = d_offsets[i+1];
2729  for (int c = 0; c < vd; ++c)
2730  {
2731  const double dofValue = d_x(t?c:i,t?i:c);
2732  for (int j = offset; j < nextOffset; ++j)
2733  {
2734  const int idx_j = d_indices[j];
2735  d_y(idx_j % nd, c, idx_j / nd) = dofValue;
2736  }
2737  }
2738  });
2739 }
2740 
2742 {
2743  // Assumes all elements have the same number of dofs
2744  const int nd = dof;
2745  const int vd = vdim;
2746  const bool t = byvdim;
2747  auto d_offsets = offsets.Read();
2748  auto d_indices = indices.Read();
2749  auto d_x = Reshape(x.Read(), nd, vd, ne);
2750  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
2751  MFEM_FORALL(i, ndofs,
2752  {
2753  const int offset = d_offsets[i];
2754  const int nextOffset = d_offsets[i + 1];
2755  for (int c = 0; c < vd; ++c)
2756  {
2757  double dofValue = 0;
2758  for (int j = offset; j < nextOffset; ++j)
2759  {
2760  const int idx_j = d_indices[j];
2761  dofValue += d_x(idx_j % nd, c, idx_j / nd);
2762  }
2763  d_y(t?c:i,t?i:c) = dofValue;
2764  }
2765  });
2766 }
2767 
2768 
2770  const IntegrationRule &ir)
2771 {
2772  fespace = &fes;
2773  qspace = NULL;
2774  IntRule = &ir;
2775  use_tensor_products = true; // not implemented yet (not used)
2776 
2777  if (fespace->GetNE() == 0) { return; }
2778  const FiniteElement *fe = fespace->GetFE(0);
2779  MFEM_VERIFY(dynamic_cast<const ScalarFiniteElement*>(fe) != NULL,
2780  "Only scalar finite elements are supported");
2781 }
2782 
2784  const QuadratureSpace &qs)
2785 {
2786  fespace = &fes;
2787  qspace = &qs;
2788  IntRule = NULL;
2789  use_tensor_products = true; // not implemented yet (not used)
2790 
2791  if (fespace->GetNE() == 0) { return; }
2792  const FiniteElement *fe = fespace->GetFE(0);
2793  MFEM_VERIFY(dynamic_cast<const ScalarFiniteElement*>(fe) != NULL,
2794  "Only scalar finite elements are supported");
2795 }
2796 
2797 template<const int T_VDIM, const int T_ND, const int T_NQ>
2799  const int NE,
2800  const int vdim,
2801  const DofToQuad &maps,
2802  const Vector &e_vec,
2803  Vector &q_val,
2804  Vector &q_der,
2805  Vector &q_det,
2806  const int eval_flags)
2807 {
2808  const int nd = maps.ndof;
2809  const int nq = maps.nqpt;
2810  const int ND = T_ND ? T_ND : nd;
2811  const int NQ = T_NQ ? T_NQ : nq;
2812  const int VDIM = T_VDIM ? T_VDIM : vdim;
2813  MFEM_VERIFY(ND <= MAX_ND2D, "");
2814  MFEM_VERIFY(NQ <= MAX_NQ2D, "");
2815  MFEM_VERIFY(VDIM == 2 || !(eval_flags & DETERMINANTS), "");
2816  auto B = Reshape(maps.B.Read(), NQ, ND);
2817  auto G = Reshape(maps.G.Read(), NQ, 2, ND);
2818  auto E = Reshape(e_vec.Read(), ND, VDIM, NE);
2819  auto val = Reshape(q_val.Write(), NQ, VDIM, NE);
2820  auto der = Reshape(q_der.Write(), NQ, VDIM, 2, NE);
2821  auto det = Reshape(q_det.Write(), NQ, NE);
2822  MFEM_FORALL(e, NE,
2823  {
2824  const int ND = T_ND ? T_ND : nd;
2825  const int NQ = T_NQ ? T_NQ : nq;
2826  const int VDIM = T_VDIM ? T_VDIM : vdim;
2827  constexpr int max_ND = T_ND ? T_ND : MAX_ND2D;
2828  constexpr int max_VDIM = T_VDIM ? T_VDIM : MAX_VDIM2D;
2829  double s_E[max_VDIM*max_ND];
2830  for (int d = 0; d < ND; d++)
2831  {
2832  for (int c = 0; c < VDIM; c++)
2833  {
2834  s_E[c+d*VDIM] = E(d,c,e);
2835  }
2836  }
2837  for (int q = 0; q < NQ; ++q)
2838  {
2839  if (eval_flags & VALUES)
2840  {
2841  double ed[max_VDIM];
2842  for (int c = 0; c < VDIM; c++) { ed[c] = 0.0; }
2843  for (int d = 0; d < ND; ++d)
2844  {
2845  const double b = B(q,d);
2846  for (int c = 0; c < VDIM; c++) { ed[c] += b*s_E[c+d*VDIM]; }
2847  }
2848  for (int c = 0; c < VDIM; c++) { val(q,c,e) = ed[c]; }
2849  }
2850  if ((eval_flags & DERIVATIVES) || (eval_flags & DETERMINANTS))
2851  {
2852  // use MAX_VDIM2D to avoid "subscript out of range" warnings
2853  double D[MAX_VDIM2D*2];
2854  for (int i = 0; i < 2*VDIM; i++) { D[i] = 0.0; }
2855  for (int d = 0; d < ND; ++d)
2856  {
2857  const double wx = G(q,0,d);
2858  const double wy = G(q,1,d);
2859  for (int c = 0; c < VDIM; c++)
2860  {
2861  double s_e = s_E[c+d*VDIM];
2862  D[c+VDIM*0] += s_e * wx;
2863  D[c+VDIM*1] += s_e * wy;
2864  }
2865  }
2866  if (eval_flags & DERIVATIVES)
2867  {
2868  for (int c = 0; c < VDIM; c++)
2869  {
2870  der(q,c,0,e) = D[c+VDIM*0];
2871  der(q,c,1,e) = D[c+VDIM*1];
2872  }
2873  }
2874  if (VDIM == 2 && (eval_flags & DETERMINANTS))
2875  {
2876  // The check (VDIM == 2) should eliminate this block when VDIM is
2877  // known at compile time and (VDIM != 2).
2878  det(q,e) = D[0]*D[3] - D[1]*D[2];
2879  }
2880  }
2881  }
2882  });
2883 }
2884 
2885 template<const int T_VDIM, const int T_ND, const int T_NQ>
2887  const int NE,
2888  const int vdim,
2889  const DofToQuad &maps,
2890  const Vector &e_vec,
2891  Vector &q_val,
2892  Vector &q_der,
2893  Vector &q_det,
2894  const int eval_flags)
2895 {
2896  const int nd = maps.ndof;
2897  const int nq = maps.nqpt;
2898  const int ND = T_ND ? T_ND : nd;
2899  const int NQ = T_NQ ? T_NQ : nq;
2900  const int VDIM = T_VDIM ? T_VDIM : vdim;
2901  MFEM_VERIFY(ND <= MAX_ND3D, "");
2902  MFEM_VERIFY(NQ <= MAX_NQ3D, "");
2903  MFEM_VERIFY(VDIM == 3 || !(eval_flags & DETERMINANTS), "");
2904  auto B = Reshape(maps.B.Read(), NQ, ND);
2905  auto G = Reshape(maps.G.Read(), NQ, 3, ND);
2906  auto E = Reshape(e_vec.Read(), ND, VDIM, NE);
2907  auto val = Reshape(q_val.Write(), NQ, VDIM, NE);
2908  auto der = Reshape(q_der.Write(), NQ, VDIM, 3, NE);
2909  auto det = Reshape(q_det.Write(), NQ, NE);
2910  MFEM_FORALL(e, NE,
2911  {
2912  const int ND = T_ND ? T_ND : nd;
2913  const int NQ = T_NQ ? T_NQ : nq;
2914  const int VDIM = T_VDIM ? T_VDIM : vdim;
2915  constexpr int max_ND = T_ND ? T_ND : MAX_ND2D;
2916  constexpr int max_VDIM = T_VDIM ? T_VDIM : MAX_VDIM2D;
2917  double s_E[max_VDIM*max_ND];
2918  for (int d = 0; d < ND; d++)
2919  {
2920  for (int c = 0; c < VDIM; c++)
2921  {
2922  s_E[c+d*VDIM] = E(d,c,e);
2923  }
2924  }
2925  for (int q = 0; q < NQ; ++q)
2926  {
2927  if (eval_flags & VALUES)
2928  {
2929  double ed[max_VDIM];
2930  for (int c = 0; c < VDIM; c++) { ed[c] = 0.0; }
2931  for (int d = 0; d < ND; ++d)
2932  {
2933  const double b = B(q,d);
2934  for (int c = 0; c < VDIM; c++) { ed[c] += b*s_E[c+d*VDIM]; }
2935  }
2936  for (int c = 0; c < VDIM; c++) { val(q,c,e) = ed[c]; }
2937  }
2938  if ((eval_flags & DERIVATIVES) || (eval_flags & DETERMINANTS))
2939  {
2940  // use MAX_VDIM3D to avoid "subscript out of range" warnings
2941  double D[MAX_VDIM3D*3];
2942  for (int i = 0; i < 3*VDIM; i++) { D[i] = 0.0; }
2943  for (int d = 0; d < ND; ++d)
2944  {
2945  const double wx = G(q,0,d);
2946  const double wy = G(q,1,d);
2947  const double wz = G(q,2,d);
2948  for (int c = 0; c < VDIM; c++)
2949  {
2950  double s_e = s_E[c+d*VDIM];
2951  D[c+VDIM*0] += s_e * wx;
2952  D[c+VDIM*1] += s_e * wy;
2953  D[c+VDIM*2] += s_e * wz;
2954  }
2955  }
2956  if (eval_flags & DERIVATIVES)
2957  {
2958  for (int c = 0; c < VDIM; c++)
2959  {
2960  der(q,c,0,e) = D[c+VDIM*0];
2961  der(q,c,1,e) = D[c+VDIM*1];
2962  der(q,c,2,e) = D[c+VDIM*2];
2963  }
2964  }
2965  if (VDIM == 3 && (eval_flags & DETERMINANTS))
2966  {
2967  // The check (VDIM == 3) should eliminate this block when VDIM is
2968  // known at compile time and (VDIM != 3).
2969  det(q,e) = D[0] * (D[4] * D[8] - D[5] * D[7]) +
2970  D[3] * (D[2] * D[7] - D[1] * D[8]) +
2971  D[6] * (D[1] * D[5] - D[2] * D[4]);
2972  }
2973  }
2974  }
2975  });
2976 }
2977 
2979  const Vector &e_vec, unsigned eval_flags,
2980  Vector &q_val, Vector &q_der, Vector &q_det) const
2981 {
2982  const int ne = fespace->GetNE();
2983  if (ne == 0) { return; }
2984  const int vdim = fespace->GetVDim();
2985  const int dim = fespace->GetMesh()->Dimension();
2986  const FiniteElement *fe = fespace->GetFE(0);
2987  const IntegrationRule *ir =
2989  const DofToQuad &maps = fe->GetDofToQuad(*ir, DofToQuad::FULL);
2990  const int nd = maps.ndof;
2991  const int nq = maps.nqpt;
2992  void (*eval_func)(
2993  const int NE,
2994  const int vdim,
2995  const DofToQuad &maps,
2996  const Vector &e_vec,
2997  Vector &q_val,
2998  Vector &q_der,
2999  Vector &q_det,
3000  const int eval_flags) = NULL;
3001  if (vdim == 1)
3002  {
3003  if (dim == 2)
3004  {
3005  switch (100*nd + nq)
3006  {
3007  // Q0
3008  case 101: eval_func = &Eval2D<1,1,1>; break;
3009  case 104: eval_func = &Eval2D<1,1,4>; break;
3010  // Q1
3011  case 404: eval_func = &Eval2D<1,4,4>; break;
3012  case 409: eval_func = &Eval2D<1,4,9>; break;
3013  // Q2
3014  case 909: eval_func = &Eval2D<1,9,9>; break;
3015  case 916: eval_func = &Eval2D<1,9,16>; break;
3016  // Q3
3017  case 1616: eval_func = &Eval2D<1,16,16>; break;
3018  case 1625: eval_func = &Eval2D<1,16,25>; break;
3019  case 1636: eval_func = &Eval2D<1,16,36>; break;
3020  // Q4
3021  case 2525: eval_func = &Eval2D<1,25,25>; break;
3022  case 2536: eval_func = &Eval2D<1,25,36>; break;
3023  case 2549: eval_func = &Eval2D<1,25,49>; break;
3024  case 2564: eval_func = &Eval2D<1,25,64>; break;
3025  }
3026  if (nq >= 100 || !eval_func)
3027  {
3028  eval_func = &Eval2D<1>;
3029  }
3030  }
3031  else if (dim == 3)
3032  {
3033  switch (1000*nd + nq)
3034  {
3035  // Q0
3036  case 1001: eval_func = &Eval3D<1,1,1>; break;
3037  case 1008: eval_func = &Eval3D<1,1,8>; break;
3038  // Q1
3039  case 8008: eval_func = &Eval3D<1,8,8>; break;
3040  case 8027: eval_func = &Eval3D<1,8,27>; break;
3041  // Q2
3042  case 27027: eval_func = &Eval3D<1,27,27>; break;
3043  case 27064: eval_func = &Eval3D<1,27,64>; break;
3044  // Q3
3045  case 64064: eval_func = &Eval3D<1,64,64>; break;
3046  case 64125: eval_func = &Eval3D<1,64,125>; break;
3047  case 64216: eval_func = &Eval3D<1,64,216>; break;
3048  // Q4
3049  case 125125: eval_func = &Eval3D<1,125,125>; break;
3050  case 125216: eval_func = &Eval3D<1,125,216>; break;
3051  }
3052  if (nq >= 1000 || !eval_func)
3053  {
3054  eval_func = &Eval3D<1>;
3055  }
3056  }
3057  }
3058  else if (vdim == dim)
3059  {
3060  if (dim == 2)
3061  {
3062  switch (100*nd + nq)
3063  {
3064  // Q1
3065  case 404: eval_func = &Eval2D<2,4,4>; break;
3066  case 409: eval_func = &Eval2D<2,4,9>; break;
3067  // Q2
3068  case 909: eval_func = &Eval2D<2,9,9>; break;
3069  case 916: eval_func = &Eval2D<2,9,16>; break;
3070  // Q3
3071  case 1616: eval_func = &Eval2D<2,16,16>; break;
3072  case 1625: eval_func = &Eval2D<2,16,25>; break;
3073  case 1636: eval_func = &Eval2D<2,16,36>; break;
3074  // Q4
3075  case 2525: eval_func = &Eval2D<2,25,25>; break;
3076  case 2536: eval_func = &Eval2D<2,25,36>; break;
3077  case 2549: eval_func = &Eval2D<2,25,49>; break;
3078  case 2564: eval_func = &Eval2D<2,25,64>; break;
3079  }
3080  if (nq >= 100 || !eval_func)
3081  {
3082  eval_func = &Eval2D<2>;
3083  }
3084  }
3085  else if (dim == 3)
3086  {
3087  switch (1000*nd + nq)
3088  {
3089  // Q1
3090  case 8008: eval_func = &Eval3D<3,8,8>; break;
3091  case 8027: eval_func = &Eval3D<3,8,27>; break;
3092  // Q2
3093  case 27027: eval_func = &Eval3D<3,27,27>; break;
3094  case 27064: eval_func = &Eval3D<3,27,64>; break;
3095  // Q3
3096  case 64064: eval_func = &Eval3D<3,64,64>; break;
3097  case 64125: eval_func = &Eval3D<3,64,125>; break;
3098  case 64216: eval_func = &Eval3D<3,64,216>; break;
3099  // Q4
3100  case 125125: eval_func = &Eval3D<3,125,125>; break;
3101  case 125216: eval_func = &Eval3D<3,125,216>; break;
3102  }
3103  if (nq >= 1000 || !eval_func)
3104  {
3105  eval_func = &Eval3D<3>;
3106  }
3107  }
3108  }
3109  if (eval_func)
3110  {
3111  eval_func(ne, vdim, maps, e_vec, q_val, q_der, q_det, eval_flags);
3112  }
3113  else
3114  {
3115  MFEM_ABORT("case not supported yet");
3116  }
3117 }
3118 
3120  unsigned eval_flags, const Vector &q_val, const Vector &q_der,
3121  Vector &e_vec) const
3122 {
3123  MFEM_ABORT("this method is not implemented yet");
3124 }
3125 
3126 } // namespace mfem
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:237
Abstract class for Finite Elements.
Definition: fe.hpp:229
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:309
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:278
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
Definition: fespace.cpp:2259
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:2715
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:359
int Size() const
Logical size of the array.
Definition: array.hpp:118
void GetFaceEdges(int i, Array< int > &, Array< int > &) const
Definition: mesh.cpp:4096
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:252
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:498
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:347
void GetVertexVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:196
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:2244
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:96
const Vector & GetWeights() const
Definition: nurbs.hpp:373
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:3853
void Load(std::istream &in, int fmt=0)
Read an Array from the stream in using format fmt. The format fmt can be:
Definition: array.cpp:53
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:668
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:344
int * GetJ()
Definition: table.hpp:114
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:105
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:85
DenseMatrix & GetPointMat()
Read and write access to the underlying point matrix describing the transformation.
Definition: eltrans.hpp:313
const Array< int > & GetMaster() const
Definition: nurbs.hpp:317
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:202
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:143
For scalar fields; preserves volume integrals.
Definition: fe.hpp:274
void GetBdrElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of bdr element i.
Definition: mesh.cpp:4064
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:769
const CoarseFineTransformations & GetDerefinementTransforms()
Definition: ncmesh.cpp:3156
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:1942
void SetSize(int dim, int connections_per_row)
Set the size and the number of connections for the table.
Definition: table.cpp:128
void LoadBE(int i, const FiniteElement *BE) const
Definition: nurbs.cpp:2779
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:877
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:4126
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
Definition: fespace.hpp:767
void BuildElementToDofTable() const
Definition: fespace.cpp:214
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:660
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:130
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
Definition: fespace.cpp:1033
void Prolongate(const Vector &x, Vector &y) const
Definition: fespace.cpp:2596
static const int NumGeom
Definition: geom.hpp:41
static const int MAX_ND2D
Definition: fespace.hpp:921
FiniteElementSpace & dom_fes
Domain FE space.
Definition: fespace.hpp:675
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:1835
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:367
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:4018
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
Definition: hypre.cpp:1086
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:400
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Definition: fespace.cpp:926
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition: fespace.cpp:1192
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:172
Ordering::Type ordering
Definition: fespace.hpp:102
const Geometry::Type geom
Definition: ex1.cpp:40
static int DecodeDof(int dof, double &sign)
Helper to remove encoded sign from a DOF.
Definition: fespace.hpp:144
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:299
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:265
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor localP[]) const
Definition: fespace.cpp:840
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:726
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:186
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:42
void GetElementFaces(int i, Array< int > &, Array< int > &) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:4258
int SizeK() const
Definition: densemat.hpp:695
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
static const int MAX_NQ3D
Definition: fespace.hpp:924
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:91
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const =0
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:53
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:795
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:472
int GetFaceOrder(int i) const
Returns the order of the i&#39;th face finite element.
Definition: fespace.cpp:101
T * GetData()
Returns the data.
Definition: array.hpp:92
int GetNKV() const
Definition: nurbs.hpp:342
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:315
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Definition: fespace.cpp:297
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:821
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe.hpp:163
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:97
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:191
OperatorHandle L2E_lex
Definition: fespace.hpp:130
void GetEntityDofs(int entity, int index, Array< int > &dofs) const
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition: fespace.cpp:571
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
Definition: fespace.hpp:681
void Save(std::ostream &out, int fmt=0) const
Save the Array to the stream out using the format fmt. The format fmt can be:
Definition: array.cpp:40
int GetBdrAttribute(int i) const
Definition: fespace.hpp:412
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe.hpp:159
void SetSize(int i, int j, int k)
Definition: densemat.hpp:699
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:366
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
void OrientedPointMatrix(DenseMatrix &oriented_matrix) const
Return the point matrix oriented according to the master and slave edges.
Definition: ncmesh.cpp:2289
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base g...
Definition: fe.cpp:111
T Sum()
Sum all entries.
Definition: array.cpp:115
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:754
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Evaluate the derivatives at quadrature points.
Definition: fespace.hpp:932
const FiniteElement * GetFaceElement(int i) const
Definition: fespace.cpp:1807
virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh.
Definition: fespace.cpp:1901
int GetMapType() const
Definition: fe.hpp:331
void GetInverseMatrix(DenseMatrix &Ainv) const
Compute and return the inverse matrix in Ainv.
Definition: densemat.cpp:4262
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller...
Definition: fespace.cpp:2088
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: fespace.cpp:2624
void Factor()
Factor the current DenseMatrix, *a.
Definition: densemat.cpp:4250
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:334
int HasFaceDofs(Geometry::Type GeomType) const
Definition: fe_coll.cpp:25
static bool IsEnabled()
Return true if any backend other than Backend::CPU is enabled.
Definition: device.hpp:195
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:96
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2467
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:313
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:159
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:94
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:127
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:759
int Size_of_connections() const
Definition: table.hpp:98
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A * B.
Definition: hypre.cpp:1543
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition: fespace.cpp:1166
bool HasGeometry(Geometry::Type geom) const
Return true iff the given geom is encountered in the mesh. Geometries of dimensions lower than Dimens...
Definition: mesh.hpp:771
void skip_comment_lines(std::istream &is, const char comment_char)
Definition: text.hpp:27
void DeleteAll()
Delete whole array.
Definition: array.hpp:785
int GetNDof() const
Definition: nurbs.hpp:352
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:756
static const int MAX_ND3D
Definition: fespace.hpp:925
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition: ncmesh.hpp:205
Array< int > dof_elem_array
Definition: fespace.hpp:113
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition: fespace.cpp:385
Array< int > indices
Definition: fespace.hpp:891
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1766
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1641
FiniteElementSpace & ran_fes
Range FE space.
Definition: fespace.hpp:676
void ConvertFromConformingVDofs(const Array< int > &cdofs, Array< int > &dofs)
Definition: fespace.cpp:422
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:789
const DenseTensor & GetPointMatrices(Geometry::Type geom) const
Definition: ncmesh.cpp:25
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:135
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication with the inverse of dense matrix.
Definition: densemat.cpp:4290
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:776
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:308
ID for class SparseMatrix.
Definition: operator.hpp:138
const Array< int > & GetSlave() const
Definition: nurbs.hpp:319
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1135
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Definition: vector.cpp:51
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:374
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:240
long GetSequence() const
Definition: mesh.hpp:1149
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1729
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: fespace.cpp:2412
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:696
General product operator: x -&gt; (A*B)(x) = A(B(x)).
Definition: operator.hpp:343
ID for the base class Operator, i.e. any type.
Definition: operator.hpp:137
int dim
Definition: ex3.cpp:48
const FiniteElement * GetEdgeElement(int i) const
Definition: fespace.cpp:1830
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Extract all column indices and values from a given row.
Definition: sparsemat.cpp:2380
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:984
void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += B * C.
Definition: densemat.cpp:3100
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:383
Data type sparse matrix.
Definition: sparsemat.hpp:40
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:690
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:272
OperatorHandle F
Forward, coarse-to-fine, operator.
Definition: fespace.hpp:770
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:146
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
Definition: fespace.hpp:909
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
Definition: fespace.cpp:2244
std::vector< Master > masters
Definition: ncmesh.hpp:189
Type
Ordering methods:
Definition: fespace.hpp:30
SparseMatrix * cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:123
void GetRow(int r, Vector &row) const
Definition: densemat.cpp:2380
const IntegrationRule * IntRule
Not owned.
Definition: fespace.hpp:916
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:463
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:431
void AddConnection(int r, int c)
Definition: table.hpp:80
void GetCoarseToFineMap(const Mesh &fine_mesh, Table &coarse_to_fine, Array< int > &coarse_to_ref_type, Table &ref_type_to_matrix, Array< Geometry::Type > &ref_type_to_geom) const
Definition: ncmesh.cpp:66
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:272
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:288
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1439
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:582
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:1341
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition: array.hpp:136
Array< Embedding > embeddings
Fine element positions in their parents.
Definition: ncmesh.hpp:64
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:261
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:57
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition: fespace.cpp:2454
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I)
Definition: fespace.cpp:532
Array< int > dof_ldof_array
Definition: fespace.hpp:113
const Array< int > & GetOrders() const
Read-only access to the orders of all knot vectors.
Definition: nurbs.hpp:337
bool Conforming() const
Definition: fespace.hpp:278
SparseMatrix * cP
Definition: fespace.hpp:121
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
Definition: fespace.cpp:2641
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:336
int Dimension() const
Definition: mesh.hpp:713
static const int MAX_VDIM2D
Definition: fespace.hpp:922
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
Definition: ncmesh.cpp:3616
Array< QuadratureInterpolator * > E2Q_array
Definition: fespace.hpp:132
void GetEdgeDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1700
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
Definition: fespace.hpp:130
int slaves_end
slave faces
Definition: ncmesh.hpp:165
void GetTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes...
Definition: fespace.cpp:1876
int SizeI() const
Definition: densemat.hpp:693
OperatorHandle B
Backward, fine-to-coarse, operator.
Definition: fespace.hpp:771
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1741
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition: fespace.cpp:902
int GetOrder() const
If all orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.
Definition: nurbs.hpp:340
List of mesh geometries stored as Array&lt;Geometry::Type&gt;.
Definition: mesh.hpp:786
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: fespace.cpp:2630
void GetBdrElementFace(int i, int *, int *) const
Return the index and the orientation of the face of bdr element i. (3D)
Definition: mesh.cpp:4296
const FiniteElementSpace & fes
Definition: fespace.hpp:883
void SetMassIntegrator(BilinearFormIntegrator *mass_integ_, bool own_mass_integ_=true)
Assign a mass integrator to be used in the construction of the backward, fine-to-coarse, transfer operator.
Definition: fespace.cpp:2370
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:980
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:398
QuadratureInterpolator(const FiniteElementSpace &fes, const IntegrationRule &ir)
Definition: fespace.cpp:2769
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:135
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const =0
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:23
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, int old_ndofs)
Definition: fespace.cpp:943
void LoadFE(int i, const FiniteElement *FE) const
Definition: nurbs.cpp:2758
void SetOrder(int Order) const
Set the order and the name, based on the given Order: either a positive number for fixed order...
Definition: fe_coll.cpp:2492
void SetIdentityTransformation(Geometry::Type GeomType)
Definition: eltrans.cpp:370
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Definition: fespace.hpp:880
std::vector< Slave > slaves
Definition: ncmesh.hpp:190
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:647
Nonconforming edge/face within a bigger edge/face.
Definition: ncmesh.hpp:172
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:311
virtual ~FiniteElementSpace()
Definition: fespace.cpp:1841
bool OwnsOperator() const
Return true if the OperatorHandle owns the held Operator.
Definition: handle.hpp:109
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:41
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Definition: fespace.cpp:494
static void Eval2D(const int NE, const int vdim, const DofToQuad &maps, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, const int eval_flags)
Template compute kernel for 2D.
Definition: fespace.cpp:2798
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:618
L2Prolongation * B
Backward, fine-to-coarse, operator.
Definition: fespace.hpp:862
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:517
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:95
void Save(std::ostream &out) const
Write the QuadratureSpace to the stream out.
Definition: fespace.cpp:2236
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:116
Table * GetElementDofTable()
Definition: nurbs.hpp:363
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe.hpp:125
void MakeVDimMatrix(SparseMatrix &mat) const
Definition: fespace.cpp:741
Table * GetBdrElementDofTable()
Definition: nurbs.hpp:364
bool own_mass_integ
Ownership flag for mass_integ.
Definition: fespace.hpp:768
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:26
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:208
void ShiftUpI()
Definition: table.cpp:119
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:7708
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:24
L2Projection * F
Forward, coarse-to-fine, operator.
Definition: fespace.hpp:861
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1143
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe.hpp:174
const Array< int > & GetDofMap() const
Get an Array&lt;int&gt; that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe.hpp:1825
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:47
bool Nonconforming() const
Definition: fespace.hpp:279
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Definition: fe.cpp:206
int SizeJ() const
Definition: densemat.hpp:694
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1558
Assuming the derivative at quadrature points form a matrix, this flag can be used to compute and stor...
Definition: fespace.hpp:936
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:181
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:673
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: fespace.cpp:2379
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
Definition: fespace.cpp:414
virtual const char * Name() const
Definition: fe_coll.hpp:53
static void Eval3D(const int NE, const int vdim, const DofToQuad &maps, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, const int eval_flags)
Template compute kernel for 3D.
Definition: fespace.cpp:2886
Class for integration point with weight.
Definition: intrules.hpp:25
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
Full multidimensional representation which does not use tensor product structure. The ordering of the...
Definition: fe.hpp:143
bool Parallel() const
Definition: fespace.hpp:689
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
Definition: mesh.cpp:4042
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:91
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:326
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:243
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:162
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:816
static const int DimStart[MaxDim+2]
Definition: geom.hpp:47
Array< int > offsets
Definition: fespace.hpp:890
int GetNConformingDofs() const
Definition: fespace.cpp:783
void ConnectBoundaries()
Definition: nurbs.cpp:1640
static const int MAX_VDIM3D
Definition: fespace.hpp:926
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1754
void Mult(const Vector &e_vec, unsigned eval_flags, Vector &q_val, Vector &q_der, Vector &q_det) const
Interpolate the E-vector e_vec to quadrature points.
Definition: fespace.cpp:2978
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition: fespace.cpp:812
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:63
General triple product operator x -&gt; A*B*C*x, with ownership of the factors.
Definition: operator.hpp:390
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Definition: fe.cpp:117
void SetFE(const FiniteElement *FE)
Definition: eltrans.hpp:299
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe.hpp:195
IsoparametricTransformation Transf
Definition: eltrans.hpp:339
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1541
Lexicographic ordering for tensor-product FiniteElements.
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
int parent
Element index in the coarse mesh.
Definition: ncmesh.hpp:48
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition: fespace.cpp:59
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: fespace.hpp:648
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition: fespace.cpp:1295
NCMesh * ncmesh
Optional non-conforming mesh extension.
Definition: mesh.hpp:182
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:764
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th face element (2D and 3D).
Definition: fespace.cpp:184
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2420
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:190
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:248
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Definition: fespace.cpp:403
virtual void Mult(const Vector &x, Vector &y) const
Perform the L2 projection onto the LOR space.
Definition: fespace.cpp:2570
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:682
void filter_dos(std::string &line)
Definition: text.hpp:41
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:803
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:685
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:3725
Vector data type.
Definition: vector.hpp:48
For scalar fields; preserves point values.
Definition: fe.hpp:273
ID for class HypreParMatrix.
Definition: operator.hpp:139
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
Definition: densemat.cpp:2649
const Table & GetElementToDofTable() const
Definition: fespace.hpp:481
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:59
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
Definition: fe.cpp:123
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition: ncmesh.hpp:212
static const int MAX_NQ2D
Definition: fespace.hpp:920
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:556
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:628
void UseExternalData(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:66
int RowSize(int i) const
Definition: table.hpp:108
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1781
NURBSExtension * NURBSext
Definition: fespace.hpp:115
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:178
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:292
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:3617
int index
Mesh number.
Definition: ncmesh.hpp:153
BiLinear2DFiniteElement QuadrilateralFE
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:656
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:368
void Save(std::ostream &out) const
Definition: fespace.cpp:2019
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
Definition: mesh.hpp:807
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:803
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:159
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:124
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:25
Integrator for (Q u, v) for VectorFiniteElements.
void MultTranspose(unsigned eval_flags, const Vector &q_val, const Vector &q_der, Vector &e_vec) const
Perform the transpose operation of Mult(). (TODO)
Definition: fespace.cpp:3119
Defines the position of a fine element within a coarse element.
Definition: ncmesh.hpp:45
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:137
int matrix
Index into the DenseTensor corresponding to the parent Geometry::Type stored in CoarseFineTransformat...
Definition: ncmesh.hpp:51
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:134
Evaluate the values at quadrature points.
Definition: fespace.hpp:931
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: fespace.cpp:2741
const QuadratureSpace * qspace
Not owned.
Definition: fespace.hpp:915
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:107