MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
gridfunc.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 GridFunction
13 
14 #include <limits>
15 #include <cstring>
16 #include <string>
17 #include <cmath>
18 #include <iostream>
19 #include <algorithm>
20 #include "fem.hpp"
21 
22 namespace mfem
23 {
24 
25 using namespace std;
26 
27 GridFunction::GridFunction(Mesh *m, std::istream &input)
28  : Vector()
29 {
30  const int bufflen = 256;
31  char buff[bufflen];
32  int vdim;
33 
34  input >> std::ws;
35  input.getline(buff, bufflen); // 'FiniteElementSpace'
36  if (strcmp(buff, "FiniteElementSpace"))
37  mfem_error("GridFunction::GridFunction():"
38  " input stream is not a GridFunction!");
39  input.getline(buff, bufflen, ' '); // 'FiniteElementCollection:'
40  input >> std::ws;
41  input.getline(buff, bufflen);
43  input.getline(buff, bufflen, ' '); // 'VDim:'
44  input >> vdim;
45  input.getline(buff, bufflen, ' '); // 'Ordering:'
46  int ordering;
47  input >> ordering;
48  input.getline(buff, bufflen); // read the empty line
49  fes = new FiniteElementSpace(m, fec, vdim, ordering);
50  Vector::Load(input, fes->GetVSize());
51 }
52 
53 GridFunction::GridFunction(Mesh *m, GridFunction *gf_array[], int num_pieces)
54 {
55  // all GridFunctions must have the same FE collection, vdim, ordering
56  int vdim, ordering;
57 
58  fes = gf_array[0]->FESpace();
60  vdim = fes->GetVDim();
61  ordering = fes->GetOrdering();
62  fes = new FiniteElementSpace(m, fec, vdim, ordering);
63  SetSize(fes->GetVSize());
64 
65  if (m->NURBSext)
66  {
67  m->NURBSext->MergeGridFunctions(gf_array, num_pieces, *this);
68  return;
69  }
70 
71  int g_ndofs = fes->GetNDofs();
72  int g_nvdofs = fes->GetNVDofs();
73  int g_nedofs = fes->GetNEDofs();
74  int g_nfdofs = fes->GetNFDofs();
75  int g_nddofs = g_ndofs - (g_nvdofs + g_nedofs + g_nfdofs);
76  int vi, ei, fi, di;
77  vi = ei = fi = di = 0;
78  for (int i = 0; i < num_pieces; i++)
79  {
80  FiniteElementSpace *l_fes = gf_array[i]->FESpace();
81  int l_ndofs = l_fes->GetNDofs();
82  int l_nvdofs = l_fes->GetNVDofs();
83  int l_nedofs = l_fes->GetNEDofs();
84  int l_nfdofs = l_fes->GetNFDofs();
85  int l_nddofs = l_ndofs - (l_nvdofs + l_nedofs + l_nfdofs);
86  const double *l_data = gf_array[i]->GetData();
87  double *g_data = data;
88  if (ordering == Ordering::byNODES)
89  {
90  for (int d = 0; d < vdim; d++)
91  {
92  memcpy(g_data+vi, l_data, l_nvdofs*sizeof(double));
93  l_data += l_nvdofs;
94  g_data += g_nvdofs;
95  memcpy(g_data+ei, l_data, l_nedofs*sizeof(double));
96  l_data += l_nedofs;
97  g_data += g_nedofs;
98  memcpy(g_data+fi, l_data, l_nfdofs*sizeof(double));
99  l_data += l_nfdofs;
100  g_data += g_nfdofs;
101  memcpy(g_data+di, l_data, l_nddofs*sizeof(double));
102  l_data += l_nddofs;
103  g_data += g_nddofs;
104  }
105  }
106  else
107  {
108  memcpy(g_data+vdim*vi, l_data, vdim*l_nvdofs*sizeof(double));
109  l_data += vdim*l_nvdofs;
110  g_data += vdim*g_nvdofs;
111  memcpy(g_data+vdim*ei, l_data, vdim*l_nedofs*sizeof(double));
112  l_data += vdim*l_nedofs;
113  g_data += vdim*g_nedofs;
114  memcpy(g_data+vdim*fi, l_data, vdim*l_nfdofs*sizeof(double));
115  l_data += vdim*l_nfdofs;
116  g_data += vdim*g_nfdofs;
117  memcpy(g_data+vdim*di, l_data, vdim*l_nddofs*sizeof(double));
118  l_data += vdim*l_nddofs;
119  g_data += vdim*g_nddofs;
120  }
121  vi += l_nvdofs;
122  ei += l_nedofs;
123  fi += l_nfdofs;
124  di += l_nddofs;
125  }
126 }
127 
129 {
130  if (fec)
131  {
132  delete fes;
133  delete fec;
134  }
135 }
136 
138 {
139  if (fec)
140  {
141  delete fes;
142  delete fec;
143  fec = NULL;
144  }
145  fes = f;
146  SetSize(fes->GetVSize());
147 }
148 
150 {
151  MFEM_ASSERT(v.Size() >= v_offset + f->GetVSize(), "");
152  if (fec)
153  {
154  delete fes;
155  delete fec;
156  fec = NULL;
157  }
158  fes = f;
159  NewDataAndSize((double *)v + v_offset, fes->GetVSize());
160 }
161 
163  GridFunction &flux,
164  Array<int>& count,
165  int wcoef,
166  int subdomain)
167 {
168  GridFunction &u = *this;
169 
170  ElementTransformation *Transf;
171 
172  FiniteElementSpace *ufes = u.FESpace();
173  FiniteElementSpace *ffes = flux.FESpace();
174 
175  int nfe = ufes->GetNE();
176  Array<int> udofs;
177  Array<int> fdofs;
178  Vector ul, fl;
179 
180  flux = 0.0;
181  count = 0;
182 
183  for (int i = 0; i < nfe; i++)
184  {
185  if (subdomain >= 0 && ufes->GetAttribute(i) != subdomain)
186  {
187  continue;
188  }
189 
190  ufes->GetElementVDofs(i, udofs);
191  ffes->GetElementVDofs(i, fdofs);
192 
193  u.GetSubVector(udofs, ul);
194 
195  Transf = ufes->GetElementTransformation(i);
196  blfi.ComputeElementFlux(*ufes->GetFE(i), *Transf, ul,
197  *ffes->GetFE(i), fl, wcoef);
198 
199  flux.AddElementVector(fdofs, fl);
200 
202  for (int j = 0; j < fdofs.Size(); j++)
203  {
204  count[fdofs[j]]++;
205  }
206  }
207 }
208 
210  GridFunction &flux, int wcoef,
211  int subdomain)
212 {
213  Array<int> count(flux.Size());
214 
215  SumFluxAndCount(blfi, flux, count, 0, subdomain);
216 
217  // complete averaging
218  for (int i = 0; i < count.Size(); i++)
219  {
220  if (count[i] != 0) { flux(i) /= count[i]; }
221  }
222 }
223 
225 {
226  if (!fes->GetNE())
227  {
228  return 0;
229  }
230 
231  const FiniteElement *fe = fes->GetFE(0);
232 
233  if (fe->GetRangeType() == FiniteElement::SCALAR)
234  {
235  return fes->GetVDim();
236  }
237  return fes->GetMesh()->SpaceDimension();
238 }
239 
240 void GridFunction::GetNodalValues(int i, Array<double> &nval, int vdim) const
241 {
242  Array<int> vdofs;
243 
244  int k;
245 
246  fes->GetElementVDofs(i, vdofs);
247  const FiniteElement *FElem = fes->GetFE(i);
248  const IntegrationRule *ElemVert =
250  int dof = FElem->GetDof();
251  int n = ElemVert->GetNPoints();
252  nval.SetSize(n);
253  vdim--;
254  Vector loc_data;
255  GetSubVector(vdofs, loc_data);
256 
257  if (FElem->GetRangeType() == FiniteElement::SCALAR)
258  {
259  Vector shape(dof);
260  for (k = 0; k < n; k++)
261  {
262  FElem->CalcShape(ElemVert->IntPoint(k), shape);
263  nval[k] = shape * ((const double *)loc_data + dof * vdim);
264  }
265  }
266  else
267  {
269  DenseMatrix vshape(dof, FElem->GetDim());
270  for (k = 0; k < n; k++)
271  {
272  Tr->SetIntPoint(&ElemVert->IntPoint(k));
273  FElem->CalcVShape(*Tr, vshape);
274  nval[k] = loc_data * (&vshape(0,vdim));
275  }
276  }
277 }
278 
279 double GridFunction::GetValue(int i, const IntegrationPoint &ip, int vdim)
280 const
281 {
282  Array<int> dofs;
283  fes->GetElementDofs(i, dofs);
284  fes->DofsToVDofs(vdim-1, dofs);
285  Vector DofVal(dofs.Size()), LocVec;
286  fes->GetFE(i)->CalcShape(ip, DofVal);
287  GetSubVector(dofs, LocVec);
288 
289  return (DofVal * LocVec);
290 }
291 
293  Vector &val) const
294 {
295  const FiniteElement *FElem = fes->GetFE(i);
296  int dof = FElem->GetDof();
297  Array<int> vdofs;
298  fes->GetElementVDofs(i, vdofs);
299  Vector loc_data;
300  GetSubVector(vdofs, loc_data);
301  if (FElem->GetRangeType() == FiniteElement::SCALAR)
302  {
303  Vector shape(dof);
304  FElem->CalcShape(ip, shape);
305  int vdim = fes->GetVDim();
306  val.SetSize(vdim);
307  for (int k = 0; k < vdim; k++)
308  {
309  val(k) = shape * ((const double *)loc_data + dof * k);
310  }
311  }
312  else
313  {
314  int spaceDim = fes->GetMesh()->SpaceDimension();
315  DenseMatrix vshape(dof, spaceDim);
317  Tr->SetIntPoint(&ip);
318  FElem->CalcVShape(*Tr, vshape);
319  val.SetSize(spaceDim);
320  vshape.MultTranspose(loc_data, val);
321  }
322 }
323 
324 void GridFunction::GetValues(int i, const IntegrationRule &ir, Vector &vals,
325  int vdim)
326 const
327 {
328  Array<int> dofs;
329  int n = ir.GetNPoints();
330  vals.SetSize(n);
331  fes->GetElementDofs(i, dofs);
332  fes->DofsToVDofs(vdim-1, dofs);
333  const FiniteElement *FElem = fes->GetFE(i);
334  int dof = FElem->GetDof();
335  Vector DofVal(dof), loc_data(dof);
336  GetSubVector(dofs, loc_data);
337  for (int k = 0; k < n; k++)
338  {
339  FElem->CalcShape(ir.IntPoint(k), DofVal);
340  vals(k) = DofVal * loc_data;
341  }
342 }
343 
344 void GridFunction::GetValues(int i, const IntegrationRule &ir, Vector &vals,
345  DenseMatrix &tr, int vdim)
346 const
347 {
349  ET = fes->GetElementTransformation(i);
350  ET->Transform(ir, tr);
351 
352  GetValues(i, ir, vals, vdim);
353 }
354 
355 int GridFunction::GetFaceValues(int i, int side, const IntegrationRule &ir,
356  Vector &vals, DenseMatrix &tr,
357  int vdim) const
358 {
359  int n, dir;
361 
362  n = ir.GetNPoints();
363  IntegrationRule eir(n); // ---
364  if (side == 2) // automatic choice of side
365  {
366  Transf = fes->GetMesh()->GetFaceElementTransformations(i, 0);
367  if (Transf->Elem2No < 0 ||
368  fes->GetAttribute(Transf->Elem1No) <=
369  fes->GetAttribute(Transf->Elem2No))
370  {
371  dir = 0;
372  }
373  else
374  {
375  dir = 1;
376  }
377  }
378  else
379  {
380  if (side == 1 && !fes->GetMesh()->FaceIsInterior(i))
381  {
382  dir = 0;
383  }
384  else
385  {
386  dir = side;
387  }
388  }
389  if (dir == 0)
390  {
391  Transf = fes->GetMesh()->GetFaceElementTransformations(i, 4);
392  Transf->Loc1.Transform(ir, eir);
393  GetValues(Transf->Elem1No, eir, vals, tr, vdim);
394  }
395  else
396  {
397  Transf = fes->GetMesh()->GetFaceElementTransformations(i, 8);
398  Transf->Loc2.Transform(ir, eir);
399  GetValues(Transf->Elem2No, eir, vals, tr, vdim);
400  }
401 
402  return dir;
403 }
404 
406  const IntegrationRule &ir,
407  DenseMatrix &vals) const
408 {
409  const FiniteElement *FElem = fes->GetFE(T.ElementNo);
410  int dof = FElem->GetDof();
411  Array<int> vdofs;
412  fes->GetElementVDofs(T.ElementNo, vdofs);
413  Vector loc_data;
414  GetSubVector(vdofs, loc_data);
415  int nip = ir.GetNPoints();
416  if (FElem->GetRangeType() == FiniteElement::SCALAR)
417  {
418  Vector shape(dof);
419  int vdim = fes->GetVDim();
420  vals.SetSize(vdim, nip);
421  for (int j = 0; j < nip; j++)
422  {
423  const IntegrationPoint &ip = ir.IntPoint(j);
424  FElem->CalcShape(ip, shape);
425  for (int k = 0; k < vdim; k++)
426  {
427  vals(k,j) = shape * ((const double *)loc_data + dof * k);
428  }
429  }
430  }
431  else
432  {
433  int spaceDim = fes->GetMesh()->SpaceDimension();
434  DenseMatrix vshape(dof, spaceDim);
435  vals.SetSize(spaceDim, nip);
436  Vector val_j;
437  for (int j = 0; j < nip; j++)
438  {
439  const IntegrationPoint &ip = ir.IntPoint(j);
440  T.SetIntPoint(&ip);
441  FElem->CalcVShape(T, vshape);
442  vals.GetColumnReference(j, val_j);
443  vshape.MultTranspose(loc_data, val_j);
444  }
445  }
446 }
447 
449  DenseMatrix &vals, DenseMatrix &tr) const
450 {
452  Tr->Transform(ir, tr);
453 
454  GetVectorValues(*Tr, ir, vals);
455 }
456 
458  int i, int side, const IntegrationRule &ir,
459  DenseMatrix &vals, DenseMatrix &tr) const
460 {
461  int n, di;
463 
464  n = ir.GetNPoints();
465  IntegrationRule eir(n); // ---
466  Transf = fes->GetMesh()->GetFaceElementTransformations(i, 0);
467  if (side == 2)
468  {
469  if (Transf->Elem2No < 0 ||
470  fes->GetAttribute(Transf->Elem1No) <=
471  fes->GetAttribute(Transf->Elem2No))
472  {
473  di = 0;
474  }
475  else
476  {
477  di = 1;
478  }
479  }
480  else
481  {
482  di = side;
483  }
484  if (di == 0)
485  {
486  Transf = fes->GetMesh()->GetFaceElementTransformations(i, 4);
487  Transf->Loc1.Transform(ir, eir);
488  GetVectorValues(Transf->Elem1No, eir, vals, tr);
489  }
490  else
491  {
492  Transf = fes->GetMesh()->GetFaceElementTransformations(i, 8);
493  Transf->Loc2.Transform(ir, eir);
494  GetVectorValues(Transf->Elem2No, eir, vals, tr);
495  }
496 
497  return di;
498 }
499 
501 {
502  // Without averaging ...
503 
504  FiniteElementSpace *orig_fes = orig_func.FESpace();
505  Array<int> vdofs, orig_vdofs;
506  Vector shape, loc_values, orig_loc_values;
507  int i, j, d, ne, dof, odof, vdim;
508 
509  ne = fes->GetNE();
510  vdim = fes->GetVDim();
511  for (i = 0; i < ne; i++)
512  {
513  fes->GetElementVDofs(i, vdofs);
514  orig_fes->GetElementVDofs(i, orig_vdofs);
515  orig_func.GetSubVector(orig_vdofs, orig_loc_values);
516  const FiniteElement *fe = fes->GetFE(i);
517  const FiniteElement *orig_fe = orig_fes->GetFE(i);
518  dof = fe->GetDof();
519  odof = orig_fe->GetDof();
520  loc_values.SetSize(dof * vdim);
521  shape.SetSize(odof);
522  const IntegrationRule &ir = fe->GetNodes();
523  for (j = 0; j < dof; j++)
524  {
525  const IntegrationPoint &ip = ir.IntPoint(j);
526  orig_fe->CalcShape(ip, shape);
527  for (d = 0; d < vdim; d++)
528  {
529  loc_values(d*dof+j) =
530  shape * ((const double *)orig_loc_values + d * odof) ;
531  }
532  }
533  SetSubVector(vdofs, loc_values);
534  }
535 }
536 
538 {
539  // Without averaging ...
540 
541  FiniteElementSpace *orig_fes = orig_func.FESpace();
542  Array<int> vdofs, orig_vdofs;
543  Vector shape, loc_values, orig_loc_values;
544  int i, j, d, nbe, dof, odof, vdim;
545 
546  nbe = fes->GetNBE();
547  vdim = fes->GetVDim();
548  for (i = 0; i < nbe; i++)
549  {
550  fes->GetBdrElementVDofs(i, vdofs);
551  orig_fes->GetBdrElementVDofs(i, orig_vdofs);
552  orig_func.GetSubVector(orig_vdofs, orig_loc_values);
553  const FiniteElement *fe = fes->GetBE(i);
554  const FiniteElement *orig_fe = orig_fes->GetBE(i);
555  dof = fe->GetDof();
556  odof = orig_fe->GetDof();
557  loc_values.SetSize(dof * vdim);
558  shape.SetSize(odof);
559  const IntegrationRule &ir = fe->GetNodes();
560  for (j = 0; j < dof; j++)
561  {
562  const IntegrationPoint &ip = ir.IntPoint(j);
563  orig_fe->CalcShape(ip, shape);
564  for (d = 0; d < vdim; d++)
565  {
566  loc_values(d*dof+j) =
567  shape * ((const double *)orig_loc_values + d * odof);
568  }
569  }
570  SetSubVector(vdofs, loc_values);
571  }
572 }
573 
575  int i, const IntegrationRule &ir, DenseMatrix &vals,
576  DenseMatrix &tr, int comp) const
577 {
578  Array<int> vdofs;
579  ElementTransformation *transf;
580 
581  int d, j, k, n, sdim, dof, ind;
582 
583  n = ir.GetNPoints();
584  fes->GetElementVDofs(i, vdofs);
585  const FiniteElement *fe = fes->GetFE(i);
586  dof = fe->GetDof();
587  sdim = fes->GetMesh()->SpaceDimension();
588  int *dofs = &vdofs[comp*dof];
589  transf = fes->GetElementTransformation(i);
590  transf->Transform(ir, tr);
591  vals.SetSize(n, sdim);
592  DenseMatrix vshape(dof, sdim);
593  double a;
594  for (k = 0; k < n; k++)
595  {
596  const IntegrationPoint &ip = ir.IntPoint(k);
597  transf->SetIntPoint(&ip);
598  fe->CalcVShape(*transf, vshape);
599  for (d = 0; d < sdim; d++)
600  {
601  a = 0.0;
602  for (j = 0; j < dof; j++)
603  if ( (ind=dofs[j]) >= 0 )
604  {
605  a += vshape(j, d) * data[ind];
606  }
607  else
608  {
609  a -= vshape(j, d) * data[-1-ind];
610  }
611  vals(k, d) = a;
612  }
613  }
614 }
615 
617 {
619  {
620  return;
621  }
622 
623  int i, j, k;
624  int vdim = fes->GetVDim();
625  int ndofs = fes->GetNDofs();
626  double *temp = new double[size];
627 
628  k = 0;
629  for (j = 0; j < ndofs; j++)
630  for (i = 0; i < vdim; i++)
631  {
632  temp[j+i*ndofs] = data[k++];
633  }
634 
635  for (i = 0; i < size; i++)
636  {
637  data[i] = temp[i];
638  }
639 
640  delete [] temp;
641 }
642 
644 {
645  int i, k;
646  Array<int> overlap(fes->GetNV());
647  Array<int> vertices;
648  DenseMatrix vals, tr;
649 
650  val.SetSize(overlap.Size());
651  overlap = 0;
652  val = 0.0;
653 
654  comp--;
655  for (i = 0; i < fes->GetNE(); i++)
656  {
657  const IntegrationRule *ir =
659  fes->GetElementVertices(i, vertices);
660  GetVectorFieldValues(i, *ir, vals, tr);
661  for (k = 0; k < ir->GetNPoints(); k++)
662  {
663  val(vertices[k]) += vals(k, comp);
664  overlap[vertices[k]]++;
665  }
666  }
667 
668  for (i = 0; i < overlap.Size(); i++)
669  {
670  val(i) /= overlap[i];
671  }
672 }
673 
675 {
676  FiniteElementSpace *new_fes = vec_field.FESpace();
677 
678  int d, i, k, ind, dof, sdim;
679  Array<int> overlap(new_fes->GetVSize());
680  Array<int> new_vdofs;
681  DenseMatrix vals, tr;
682 
683  sdim = fes->GetMesh()->SpaceDimension();
684  overlap = 0;
685  vec_field = 0.0;
686 
687  for (i = 0; i < new_fes->GetNE(); i++)
688  {
689  const FiniteElement *fe = new_fes->GetFE(i);
690  const IntegrationRule &ir = fe->GetNodes();
691  GetVectorFieldValues(i, ir, vals, tr, comp);
692  new_fes->GetElementVDofs(i, new_vdofs);
693  dof = fe->GetDof();
694  for (d = 0; d < sdim; d++)
695  {
696  for (k = 0; k < dof; k++)
697  {
698  if ( (ind=new_vdofs[dof*d+k]) < 0 )
699  {
700  ind = -1-ind, vals(k, d) = - vals(k, d);
701  }
702  vec_field(ind) += vals(k, d);
703  overlap[ind]++;
704  }
705  }
706  }
707 
708  for (i = 0; i < overlap.Size(); i++)
709  {
710  vec_field(i) /= overlap[i];
711  }
712 }
713 
714 void GridFunction::GetDerivative(int comp, int der_comp, GridFunction &der)
715 {
716  FiniteElementSpace * der_fes = der.FESpace();
717  ElementTransformation * transf;
718  Array<int> overlap(der_fes->GetVSize());
719  Array<int> der_dofs, vdofs;
720  DenseMatrix dshape, inv_jac;
721  Vector pt_grad, loc_func;
722  int i, j, k, dim, dof, der_dof, ind;
723  double a;
724 
725  for (i = 0; i < overlap.Size(); i++)
726  {
727  overlap[i] = 0;
728  }
729  der = 0.0;
730 
731  comp--;
732  for (i = 0; i < der_fes->GetNE(); i++)
733  {
734  const FiniteElement *der_fe = der_fes->GetFE(i);
735  const FiniteElement *fe = fes->GetFE(i);
736  const IntegrationRule &ir = der_fe->GetNodes();
737  der_fes->GetElementDofs(i, der_dofs);
738  fes->GetElementVDofs(i, vdofs);
739  dim = fe->GetDim();
740  dof = fe->GetDof();
741  der_dof = der_fe->GetDof();
742  dshape.SetSize(dof, dim);
743  inv_jac.SetSize(dim);
744  pt_grad.SetSize(dim);
745  loc_func.SetSize(dof);
746  transf = fes->GetElementTransformation(i);
747  for (j = 0; j < dof; j++)
748  loc_func(j) = ( (ind=vdofs[comp*dof+j]) >= 0 ) ?
749  (data[ind]) : (-data[-1-ind]);
750  for (k = 0; k < der_dof; k++)
751  {
752  const IntegrationPoint &ip = ir.IntPoint(k);
753  fe->CalcDShape(ip, dshape);
754  dshape.MultTranspose(loc_func, pt_grad);
755  transf->SetIntPoint(&ip);
756  CalcInverse(transf->Jacobian(), inv_jac);
757  a = 0.0;
758  for (j = 0; j < dim; j++)
759  {
760  a += inv_jac(j, der_comp) * pt_grad(j);
761  }
762  der(der_dofs[k]) += a;
763  overlap[der_dofs[k]]++;
764  }
765  }
766 
767  for (i = 0; i < overlap.Size(); i++)
768  {
769  der(i) /= overlap[i];
770  }
771 }
772 
773 
776 {
777  int elNo = T.ElementNo;
778  const FiniteElement *FElem = fes->GetFE(elNo);
779  int dim = FElem->GetDim(), dof = FElem->GetDof();
780  Array<int> vdofs;
781  fes->GetElementVDofs(elNo, vdofs);
782  Vector loc_data;
783  GetSubVector(vdofs, loc_data);
784  // assuming scalar FE
785  int vdim = fes->GetVDim();
786  DenseMatrix dshape(dof, dim);
787  FElem->CalcDShape(T.GetIntPoint(), dshape);
788  gh.SetSize(vdim, dim);
789  DenseMatrix loc_data_mat(loc_data.StealData(), dof, vdim);
790  MultAtB(loc_data_mat, dshape, gh);
791 }
792 
794 {
795  double div_v;
796  int elNo = tr.ElementNo;
797  const FiniteElement *FElem = fes->GetFE(elNo);
798  if (FElem->GetRangeType() == FiniteElement::SCALAR)
799  {
800  DenseMatrix grad_hat;
801  GetVectorGradientHat(tr, grad_hat);
802  const DenseMatrix &J = tr.Jacobian();
803  DenseMatrix Jinv(J.Width(), J.Height());
804  CalcInverse(J, Jinv);
805  div_v = 0.0;
806  for (int i = 0; i < Jinv.Width(); i++)
807  for (int j = 0; j < Jinv.Height(); j++)
808  {
809  div_v += grad_hat(i, j) * Jinv(j, i);
810  }
811  }
812  else
813  {
814  // Assuming RT-type space
815  Array<int> dofs;
816  fes->GetElementDofs(elNo, dofs);
817  Vector loc_data, divshape(FElem->GetDof());
818  GetSubVector(dofs, loc_data);
819  FElem->CalcDivShape(tr.GetIntPoint(), divshape);
820  div_v = (loc_data * divshape) / tr.Weight();
821  }
822  return div_v;
823 }
824 
826 {
827  int elNo = tr.ElementNo;
828  const FiniteElement *FElem = fes->GetFE(elNo);
829  if (FElem->GetRangeType() == FiniteElement::SCALAR)
830  {
831  DenseMatrix grad_hat;
832  GetVectorGradientHat(tr, grad_hat);
833  const DenseMatrix &J = tr.Jacobian();
834  DenseMatrix Jinv(J.Width(), J.Height());
835  CalcInverse(J, Jinv);
836  DenseMatrix grad(grad_hat.Height(), Jinv.Width()); // vdim x FElem->Dim
837  Mult(grad_hat, Jinv, grad);
838  MFEM_ASSERT(grad.Height() == grad.Width(), "");
839  if (grad.Height() == 3)
840  {
841  curl.SetSize(3);
842  curl(0) = grad(2,1) - grad(1,2);
843  curl(1) = grad(0,2) - grad(2,0);
844  curl(2) = grad(1,0) - grad(0,1);
845  }
846  else if (grad.Height() == 2)
847  {
848  curl.SetSize(1);
849  curl(0) = grad(1,0) - grad(0,1);
850  }
851  }
852  else
853  {
854  // Assuming ND-type space
855  Array<int> dofs;
856  fes->GetElementDofs(elNo, dofs);
857  Vector loc_data;
858  GetSubVector(dofs, loc_data);
859  DenseMatrix curl_shape(FElem->GetDof(), FElem->GetDim() == 3 ? 3 : 1);
860  FElem->CalcCurlShape(tr.GetIntPoint(), curl_shape);
861  curl.SetSize(curl_shape.Width());
862  if (curl_shape.Width() == 3)
863  {
864  double curl_hat[3];
865  curl_shape.MultTranspose(loc_data, curl_hat);
866  tr.Jacobian().Mult(curl_hat, curl);
867  }
868  else
869  {
870  curl_shape.MultTranspose(loc_data, curl);
871  }
872  curl /= tr.Weight();
873  }
874 }
875 
877 {
878  int elNo = tr.ElementNo;
879  const FiniteElement *fe = fes->GetFE(elNo);
880  int dim = fe->GetDim(), dof = fe->GetDof();
881  DenseMatrix dshape(dof, dim), Jinv(dim);
882  Vector lval, gh(dim);
883  Array<int> dofs;
884 
885  grad.SetSize(dim);
886  fes->GetElementDofs(elNo, dofs);
887  GetSubVector(dofs, lval);
888  fe->CalcDShape(tr.GetIntPoint(), dshape);
889  dshape.MultTranspose(lval, gh);
890  CalcInverse(tr.Jacobian(), Jinv);
891  Jinv.MultTranspose(gh, grad);
892 }
893 
894 void GridFunction::GetGradients(const int elem, const IntegrationRule &ir,
895  DenseMatrix &grad)
896 {
897  const FiniteElement *fe = fes->GetFE(elem);
899  DenseMatrix dshape(fe->GetDof(), fe->GetDim());
900  DenseMatrix Jinv(fe->GetDim());
901  Vector lval, gh(fe->GetDim()), gcol;
902  Array<int> dofs;
903  fes->GetElementDofs(elem, dofs);
904  GetSubVector(dofs, lval);
905  grad.SetSize(fe->GetDim(), ir.GetNPoints());
906  for (int i = 0; i < ir.GetNPoints(); i++)
907  {
908  const IntegrationPoint &ip = ir.IntPoint(i);
909  fe->CalcDShape(ip, dshape);
910  dshape.MultTranspose(lval, gh);
911  Tr->SetIntPoint(&ip);
912  grad.GetColumnReference(i, gcol);
913  CalcInverse(Tr->Jacobian(), Jinv);
914  Jinv.MultTranspose(gh, gcol);
915  }
916 }
917 
920 {
921  DenseMatrix grad_hat;
922  GetVectorGradientHat(tr, grad_hat);
923  const DenseMatrix &J = tr.Jacobian();
924  DenseMatrix Jinv(J.Width(), J.Height());
925  CalcInverse(J, Jinv);
926  grad.SetSize(grad_hat.Height(), Jinv.Width());
927  Mult(grad_hat, Jinv, grad);
928 }
929 
931 {
932  MassIntegrator Mi;
933  DenseMatrix loc_mass;
934  Array<int> te_dofs, tr_dofs;
935  Vector loc_avgs, loc_this;
936  Vector int_psi(avgs.Size());
937 
938  avgs = 0.0;
939  int_psi = 0.0;
940  for (int i = 0; i < fes->GetNE(); i++)
941  {
942  Mi.AssembleElementMatrix2(*fes->GetFE(i), *avgs.FESpace()->GetFE(i),
943  *fes->GetElementTransformation(i), loc_mass);
944  fes->GetElementDofs(i, tr_dofs);
945  avgs.FESpace()->GetElementDofs(i, te_dofs);
946  GetSubVector(tr_dofs, loc_this);
947  loc_avgs.SetSize(te_dofs.Size());
948  loc_mass.Mult(loc_this, loc_avgs);
949  avgs.AddElementVector(te_dofs, loc_avgs);
950  loc_this = 1.0; // assume the local basis for 'this' sums to 1
951  loc_mass.Mult(loc_this, loc_avgs);
952  int_psi.AddElementVector(te_dofs, loc_avgs);
953  }
954  for (int i = 0; i < avgs.Size(); i++)
955  {
956  avgs(i) /= int_psi(i);
957  }
958 }
959 
961 {
962  // Assuming that the projection matrix is the same for all elements
963  Mesh *mesh = fes->GetMesh();
964  DenseMatrix P;
965 
966  if (!fes->GetNE())
967  {
968  return;
969  }
970 
971  fes->GetFE(0)->Project(*src.fes->GetFE(0),
972  *mesh->GetElementTransformation(0), P);
973  int vdim = fes->GetVDim();
974  if (vdim != src.fes->GetVDim())
975  mfem_error("GridFunction::ProjectGridFunction() :"
976  " incompatible vector dimensions!");
977  Array<int> src_vdofs, dest_vdofs;
978  Vector src_lvec, dest_lvec(vdim*P.Height());
979 
980  for (int i = 0; i < mesh->GetNE(); i++)
981  {
982  src.fes->GetElementVDofs(i, src_vdofs);
983  src.GetSubVector(src_vdofs, src_lvec);
984  for (int vd = 0; vd < vdim; vd++)
985  {
986  P.Mult(&src_lvec[vd*P.Width()], &dest_lvec[vd*P.Height()]);
987  }
988  fes->GetElementVDofs(i, dest_vdofs);
989  SetSubVector(dest_vdofs, dest_lvec);
990  }
991 }
992 
993 void GridFunction::ImposeBounds(int i, const Vector &weights,
994  const Vector &_lo, const Vector &_hi)
995 {
996  Array<int> vdofs;
997  fes->GetElementVDofs(i, vdofs);
998  int size = vdofs.Size();
999  Vector vals, new_vals(size);
1000  GetSubVector(vdofs, vals);
1001 
1002  MFEM_ASSERT(weights.Size() == size, "Different # of weights and dofs.");
1003  MFEM_ASSERT(_lo.Size() == size, "Different # of lower bounds and dofs.");
1004  MFEM_ASSERT(_hi.Size() == size, "Different # of upper bounds and dofs.");
1005 
1006  int max_iter = 30;
1007  double tol = 1.e-12;
1008  SLBQPOptimizer slbqp;
1009  slbqp.SetMaxIter(max_iter);
1010  slbqp.SetAbsTol(1.0e-18);
1011  slbqp.SetRelTol(tol);
1012  slbqp.SetBounds(_lo, _hi);
1013  slbqp.SetLinearConstraint(weights, weights * vals);
1014  slbqp.SetPrintLevel(0); // print messages only if not converged
1015  slbqp.Mult(vals, new_vals);
1016 
1017  SetSubVector(vdofs, new_vals);
1018 }
1019 
1020 void GridFunction::ImposeBounds(int i, const Vector &weights,
1021  double _min, double _max)
1022 {
1023  Array<int> vdofs;
1024  fes->GetElementVDofs(i, vdofs);
1025  int size = vdofs.Size();
1026  Vector vals, new_vals(size);
1027  GetSubVector(vdofs, vals);
1028 
1029  double max_val = vals.Max();
1030  double min_val = vals.Min();
1031 
1032  if (max_val <= _min)
1033  {
1034  new_vals = _min;
1035  SetSubVector(vdofs, new_vals);
1036  return;
1037  }
1038 
1039  if (_min <= min_val && max_val <= _max)
1040  {
1041  return;
1042  }
1043 
1044  Vector minv(size), maxv(size);
1045  minv = (_min > min_val) ? _min : min_val;
1046  maxv = (_max < max_val) ? _max : max_val;
1047 
1048  ImposeBounds(i, weights, minv, maxv);
1049 }
1050 
1051 void GridFunction::GetNodalValues(Vector &nval, int vdim) const
1052 {
1053  int i, j;
1054  Array<int> vertices;
1055  Array<double> values;
1056  Array<int> overlap(fes->GetNV());
1057  nval.SetSize(fes->GetNV());
1058 
1059  nval = 0.0;
1060  overlap = 0;
1061  for (i = 0; i < fes->GetNE(); i++)
1062  {
1063  fes->GetElementVertices(i, vertices);
1064  GetNodalValues(i, values, vdim);
1065  for (j = 0; j < vertices.Size(); j++)
1066  {
1067  nval(vertices[j]) += values[j];
1068  overlap[vertices[j]]++;
1069  }
1070  }
1071  for (i = 0; i < overlap.Size(); i++)
1072  {
1073  nval(i) /= overlap[i];
1074  }
1075 }
1076 
1078  double &integral)
1079 {
1080  if (!fes->GetNE())
1081  {
1082  integral = 0.0;
1083  return;
1084  }
1085 
1086  Mesh *mesh = fes->GetMesh();
1087  const int dim = mesh->Dimension();
1088  const double *center = delta_coeff.Center();
1089  const double *vert = mesh->GetVertex(0);
1090  double min_dist, dist;
1091  int v_idx = 0;
1092 
1093  // find the vertex closest to the center of the delta function
1094  min_dist = Distance(center, vert, dim);
1095  for (int i = 0; i < mesh->GetNV(); i++)
1096  {
1097  vert = mesh->GetVertex(i);
1098  dist = Distance(center, vert, dim);
1099  if (dist < min_dist)
1100  {
1101  min_dist = dist;
1102  v_idx = i;
1103  }
1104  }
1105 
1106  (*this) = 0.0;
1107  integral = 0.0;
1108 
1109  if (min_dist >= delta_coeff.Tol())
1110  {
1111  return;
1112  }
1113 
1114  // find the elements that have 'v_idx' as a vertex
1115  MassIntegrator Mi(*delta_coeff.Weight());
1116  DenseMatrix loc_mass;
1117  Array<int> vdofs, vertices;
1118  Vector vals, loc_mass_vals;
1119  for (int i = 0; i < mesh->GetNE(); i++)
1120  {
1121  mesh->GetElementVertices(i, vertices);
1122  for (int j = 0; j < vertices.Size(); j++)
1123  if (vertices[j] == v_idx)
1124  {
1125  const FiniteElement *fe = fes->GetFE(i);
1126  Mi.AssembleElementMatrix(*fe, *fes->GetElementTransformation(i),
1127  loc_mass);
1128  vals.SetSize(fe->GetDof());
1129  fe->ProjectDelta(j, vals);
1130  fes->GetElementVDofs(i, vdofs);
1131  SetSubVector(vdofs, vals);
1132  loc_mass_vals.SetSize(vals.Size());
1133  loc_mass.Mult(vals, loc_mass_vals);
1134  integral += loc_mass_vals.Sum(); // partition of unity basis
1135  break;
1136  }
1137  }
1138 }
1139 
1141 {
1142  DeltaCoefficient *delta_c = dynamic_cast<DeltaCoefficient *>(&coeff);
1143 
1144  if (delta_c == NULL)
1145  {
1146  Array<int> vdofs;
1147  Vector vals;
1148 
1149  for (int i = 0; i < fes->GetNE(); i++)
1150  {
1151  fes->GetElementVDofs(i, vdofs);
1152  vals.SetSize(vdofs.Size());
1153  fes->GetFE(i)->Project(coeff, *fes->GetElementTransformation(i), vals);
1154  SetSubVector(vdofs, vals);
1155  }
1156  }
1157  else
1158  {
1159  double integral;
1160 
1161  ProjectDeltaCoefficient(*delta_c, integral);
1162 
1163  (*this) *= (delta_c->Scale() / integral);
1164  }
1165 }
1166 
1168  Coefficient &coeff, Array<int> &dofs, int vd)
1169 {
1170  int el = -1;
1171  ElementTransformation *T = NULL;
1172  const FiniteElement *fe = NULL;
1173 
1174  for (int i = 0; i < dofs.Size(); i++)
1175  {
1176  int dof = dofs[i], j = fes->GetElementForDof(dof);
1177  if (el != j)
1178  {
1179  el = j;
1180  T = fes->GetElementTransformation(el);
1181  fe = fes->GetFE(el);
1182  }
1183  int vdof = fes->DofToVDof(dof, vd);
1184  int ld = fes->GetLocalDofForDof(dof);
1185  const IntegrationPoint &ip = fe->GetNodes().IntPoint(ld);
1186  T->SetIntPoint(&ip);
1187  (*this)(vdof) = coeff.Eval(*T, ip);
1188  }
1189 }
1190 
1192 {
1193  int i;
1194  Array<int> vdofs;
1195  Vector vals;
1196 
1197  for (i = 0; i < fes->GetNE(); i++)
1198  {
1199  fes->GetElementVDofs(i, vdofs);
1200  vals.SetSize(vdofs.Size());
1201  fes->GetFE(i)->Project(vcoeff, *fes->GetElementTransformation(i), vals);
1202  SetSubVector(vdofs, vals);
1203  }
1204 }
1205 
1207  VectorCoefficient &vcoeff, Array<int> &dofs)
1208 {
1209  int el = -1;
1210  ElementTransformation *T = NULL;
1211  const FiniteElement *fe = NULL;
1212 
1213  Vector val;
1214 
1215  for (int i = 0; i < dofs.Size(); i++)
1216  {
1217  int dof = dofs[i], j = fes->GetElementForDof(dof);
1218  if (el != j)
1219  {
1220  el = j;
1221  T = fes->GetElementTransformation(el);
1222  fe = fes->GetFE(el);
1223  }
1224  int ld = fes->GetLocalDofForDof(dof);
1225  const IntegrationPoint &ip = fe->GetNodes().IntPoint(ld);
1226  T->SetIntPoint(&ip);
1227  vcoeff.Eval(val, *T, ip);
1228  for (int vd = 0; vd < fes->GetVDim(); vd ++)
1229  {
1230  int vdof = fes->DofToVDof(dof, vd);
1231  (*this)(vdof) = val(vd);
1232  }
1233  }
1234 }
1235 
1237 {
1238  int i, j, fdof, d, ind, vdim;
1239  double val;
1240  const FiniteElement *fe;
1241  ElementTransformation *transf;
1242  Array<int> vdofs;
1243 
1244  vdim = fes->GetVDim();
1245  for (i = 0; i < fes->GetNE(); i++)
1246  {
1247  fe = fes->GetFE(i);
1248  fdof = fe->GetDof();
1249  transf = fes->GetElementTransformation(i);
1250  const IntegrationRule &ir = fe->GetNodes();
1251  fes->GetElementVDofs(i, vdofs);
1252  for (j = 0; j < fdof; j++)
1253  {
1254  const IntegrationPoint &ip = ir.IntPoint(j);
1255  transf->SetIntPoint(&ip);
1256  for (d = 0; d < vdim; d++)
1257  {
1258  val = coeff[d]->Eval(*transf, ip);
1259  if ( (ind = vdofs[fdof*d+j]) < 0 )
1260  {
1261  val = -val, ind = -1-ind;
1262  }
1263  (*this)(ind) = val;
1264  }
1265  }
1266  }
1267 }
1268 
1270  Array<int> &dof_attr)
1271 {
1272  Array<int> vdofs;
1273  Vector vals;
1274 
1275  // maximal element attribute for each dof
1276  dof_attr.SetSize(fes->GetVSize());
1277  dof_attr = -1;
1278 
1279  // local projection
1280  for (int i = 0; i < fes->GetNE(); i++)
1281  {
1282  fes->GetElementVDofs(i, vdofs);
1283  vals.SetSize(vdofs.Size());
1284  fes->GetFE(i)->Project(coeff, *fes->GetElementTransformation(i), vals);
1285 
1286  // the values in shared dofs are determined from the element with maximal
1287  // attribute
1288  int attr = fes->GetAttribute(i);
1289  for (int j = 0; j < vdofs.Size(); j++)
1290  {
1291  if (attr > dof_attr[vdofs[j]])
1292  {
1293  (*this)(vdofs[j]) = vals[j];
1294  dof_attr[vdofs[j]] = attr;
1295  }
1296  }
1297  }
1298 }
1299 
1301 {
1302  Array<int> dof_attr;
1303  ProjectDiscCoefficient(coeff, dof_attr);
1304 }
1305 
1307  Coefficient *coeff[], Array<int> &attr)
1308 {
1309  int i, j, fdof, d, ind, vdim;
1310  double val;
1311  const FiniteElement *fe;
1312  ElementTransformation *transf;
1313  Array<int> vdofs;
1314 
1315  vdim = fes->GetVDim();
1316  for (i = 0; i < fes->GetNBE(); i++)
1317  {
1318  if (attr[fes->GetBdrAttribute(i) - 1])
1319  {
1320  fe = fes->GetBE(i);
1321  fdof = fe->GetDof();
1322  transf = fes->GetBdrElementTransformation(i);
1323  const IntegrationRule &ir = fe->GetNodes();
1324  fes->GetBdrElementVDofs(i, vdofs);
1325 
1326  for (j = 0; j < fdof; j++)
1327  {
1328  const IntegrationPoint &ip = ir.IntPoint(j);
1329  transf->SetIntPoint(&ip);
1330  for (d = 0; d < vdim; d++)
1331  {
1332  val = coeff[d]->Eval(*transf, ip);
1333  if ( (ind = vdofs[fdof*d+j]) < 0 )
1334  {
1335  val = -val, ind = -1-ind;
1336  }
1337  (*this)(ind) = val;
1338  }
1339  }
1340  }
1341  }
1342 
1343  // In the case of partially conforming space, i.e. (fes->cP != NULL), we need
1344  // to set the values of all dofs on which the dofs set above depend.
1345  // Dependency is defined from the matrix A = cP.cR: dof i depends on dof j
1346  // iff A_ij != 0. It is sufficient to resolve just the first level of
1347  // dependency since A is a projection matrix: A^n = A due to cR.cP = I.
1348  // Cases like this arise in 3D when boundary edges are constraint by (depend
1349  // on) internal faces/elements.
1350  // We use the virtual method GetBoundaryClosure from NCMesh to resolve the
1351  // dependencies.
1352 
1353  if (fes->Nonconforming() && fes->GetMesh()->Dimension() == 3)
1354  {
1355  Vector vals;
1356  Mesh *mesh = fes->GetMesh();
1357  NCMesh *ncmesh = mesh->ncmesh;
1358  Array<int> bdr_edges, bdr_vertices;
1359  ncmesh->GetBoundaryClosure(attr, bdr_vertices, bdr_edges);
1360 
1361  for (i = 0; i < bdr_edges.Size(); i++)
1362  {
1363  int edge = bdr_edges[i];
1364  fes->GetEdgeVDofs(edge, vdofs);
1365  if (vdofs.Size() == 0) { continue; }
1366 
1367  transf = mesh->GetEdgeTransformation(edge);
1368  transf->Attribute = -1; // FIXME: set the boundary attribute
1369  fe = fes->GetEdgeElement(edge);
1370  vals.SetSize(fe->GetDof());
1371  for (d = 0; d < vdim; d++)
1372  {
1373  fe->Project(*coeff[d], *transf, vals);
1374  for (int k = 0; k < vals.Size(); k++)
1375  {
1376  (*this)(vdofs[d*vals.Size()+k]) = vals(k);
1377  }
1378  }
1379  }
1380  }
1381 }
1382 
1384  VectorCoefficient &vcoeff, Array<int> &bdr_attr)
1385 {
1386 #if 0
1387  // implementation for the case when the face dofs are integrals of the
1388  // normal component.
1389  const FiniteElement *fe;
1391  Array<int> dofs;
1392  int dim = vcoeff.GetVDim();
1393  Vector vc(dim), nor(dim), lvec, shape;
1394 
1395  for (int i = 0; i < fes->GetNBE(); i++)
1396  {
1397  if (bdr_attr[fes->GetBdrAttribute(i)-1] == 0)
1398  {
1399  continue;
1400  }
1401  fe = fes->GetBE(i);
1403  int intorder = 2*fe->GetOrder(); // !!!
1404  const IntegrationRule &ir = IntRules.Get(fe->GetGeomType(), intorder);
1405  int nd = fe->GetDof();
1406  lvec.SetSize(nd);
1407  shape.SetSize(nd);
1408  lvec = 0.0;
1409  for (int j = 0; j < ir.GetNPoints(); j++)
1410  {
1411  const IntegrationPoint &ip = ir.IntPoint(j);
1412  T->SetIntPoint(&ip);
1413  vcoeff.Eval(vc, *T, ip);
1414  CalcOrtho(T->Jacobian(), nor);
1415  fe->CalcShape(ip, shape);
1416  lvec.Add(ip.weight * (vc * nor), shape);
1417  }
1418  fes->GetBdrElementDofs(i, dofs);
1419  SetSubVector(dofs, lvec);
1420  }
1421 #else
1422  // implementation for the case when the face dofs are scaled point
1423  // values of the normal component.
1424  const FiniteElement *fe;
1426  Array<int> dofs;
1427  int dim = vcoeff.GetVDim();
1428  Vector vc(dim), nor(dim), lvec;
1429 
1430  for (int i = 0; i < fes->GetNBE(); i++)
1431  {
1432  if (bdr_attr[fes->GetBdrAttribute(i)-1] == 0)
1433  {
1434  continue;
1435  }
1436  fe = fes->GetBE(i);
1438  const IntegrationRule &ir = fe->GetNodes();
1439  lvec.SetSize(fe->GetDof());
1440  for (int j = 0; j < ir.GetNPoints(); j++)
1441  {
1442  const IntegrationPoint &ip = ir.IntPoint(j);
1443  T->SetIntPoint(&ip);
1444  vcoeff.Eval(vc, *T, ip);
1445  CalcOrtho(T->Jacobian(), nor);
1446  lvec(j) = (vc * nor);
1447  }
1448  fes->GetBdrElementDofs(i, dofs);
1449  SetSubVector(dofs, lvec);
1450  }
1451 #endif
1452 }
1453 
1455  VectorCoefficient &vcoeff, Array<int> &bdr_attr)
1456 {
1457  const FiniteElement *fe;
1459  Array<int> dofs;
1460  Vector lvec;
1461 
1462  for (int i = 0; i < fes->GetNBE(); i++)
1463  {
1464  if (bdr_attr[fes->GetBdrAttribute(i)-1] == 0)
1465  {
1466  continue;
1467  }
1468  fe = fes->GetBE(i);
1470  fes->GetBdrElementDofs(i, dofs);
1471  lvec.SetSize(fe->GetDof());
1472  fe->Project(vcoeff, *T, lvec);
1473  SetSubVector(dofs, lvec);
1474  }
1475 
1476  if (fes->Nonconforming() && fes->GetMesh()->Dimension() == 3)
1477  {
1478  Mesh *mesh = fes->GetMesh();
1479  NCMesh *ncmesh = mesh->ncmesh;
1480  Array<int> bdr_edges, bdr_vertices;
1481  ncmesh->GetBoundaryClosure(bdr_attr, bdr_vertices, bdr_edges);
1482 
1483  for (int i = 0; i < bdr_edges.Size(); i++)
1484  {
1485  int edge = bdr_edges[i];
1486  fes->GetEdgeDofs(edge, dofs);
1487  if (dofs.Size() == 0) { continue; }
1488 
1489  T = mesh->GetEdgeTransformation(edge);
1490  T->Attribute = -1; // FIXME: set the boundary attribute
1491  fe = fes->GetEdgeElement(edge);
1492  lvec.SetSize(fe->GetDof());
1493  fe->Project(vcoeff, *T, lvec);
1494  SetSubVector(dofs, lvec);
1495  }
1496  }
1497 }
1498 
1500  Coefficient *exsol[], const IntegrationRule *irs[]) const
1501 {
1502  double error = 0.0, a;
1503  const FiniteElement *fe;
1504  ElementTransformation *transf;
1505  Vector shape;
1506  Array<int> vdofs;
1507  int fdof, d, i, intorder, j, k;
1508 
1509  for (i = 0; i < fes->GetNE(); i++)
1510  {
1511  fe = fes->GetFE(i);
1512  fdof = fe->GetDof();
1513  transf = fes->GetElementTransformation(i);
1514  shape.SetSize(fdof);
1515  intorder = 2*fe->GetOrder() + 1; // <----------
1516  const IntegrationRule *ir;
1517  if (irs)
1518  {
1519  ir = irs[fe->GetGeomType()];
1520  }
1521  else
1522  {
1523  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
1524  }
1525  fes->GetElementVDofs(i, vdofs);
1526  for (j = 0; j < ir->GetNPoints(); j++)
1527  {
1528  const IntegrationPoint &ip = ir->IntPoint(j);
1529  fe->CalcShape(ip, shape);
1530  for (d = 0; d < fes->GetVDim(); d++)
1531  {
1532  a = 0;
1533  for (k = 0; k < fdof; k++)
1534  if (vdofs[fdof*d+k] >= 0)
1535  {
1536  a += (*this)(vdofs[fdof*d+k]) * shape(k);
1537  }
1538  else
1539  {
1540  a -= (*this)(-1-vdofs[fdof*d+k]) * shape(k);
1541  }
1542  transf->SetIntPoint(&ip);
1543  a -= exsol[d]->Eval(*transf, ip);
1544  error += ip.weight * transf->Weight() * a * a;
1545  }
1546  }
1547  }
1548 
1549  if (error < 0.0)
1550  {
1551  return -sqrt(-error);
1552  }
1553  return sqrt(error);
1554 }
1555 
1557  VectorCoefficient &exsol, const IntegrationRule *irs[],
1558  Array<int> *elems) const
1559 {
1560  double error = 0.0;
1561  const FiniteElement *fe;
1563  DenseMatrix vals, exact_vals;
1564  Vector loc_errs;
1565 
1566  for (int i = 0; i < fes->GetNE(); i++)
1567  {
1568  if (elems != NULL && (*elems)[i] == 0) { continue; }
1569  fe = fes->GetFE(i);
1570  int intorder = 2*fe->GetOrder() + 1; // <----------
1571  const IntegrationRule *ir;
1572  if (irs)
1573  {
1574  ir = irs[fe->GetGeomType()];
1575  }
1576  else
1577  {
1578  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
1579  }
1580  T = fes->GetElementTransformation(i);
1581  GetVectorValues(*T, *ir, vals);
1582  exsol.Eval(exact_vals, *T, *ir);
1583  vals -= exact_vals;
1584  loc_errs.SetSize(vals.Width());
1585  vals.Norm2(loc_errs);
1586  for (int j = 0; j < ir->GetNPoints(); j++)
1587  {
1588  const IntegrationPoint &ip = ir->IntPoint(j);
1589  T->SetIntPoint(&ip);
1590  error += ip.weight * T->Weight() * (loc_errs(j) * loc_errs(j));
1591  }
1592  }
1593 
1594  if (error < 0.0)
1595  {
1596  return -sqrt(-error);
1597  }
1598  return sqrt(error);
1599 }
1600 
1602  Coefficient *exsol, VectorCoefficient *exgrad,
1603  Coefficient *ell_coeff, double Nu, int norm_type) const
1604 {
1605  // assuming vdim is 1
1606  int i, fdof, dim, intorder, j, k;
1607  Mesh *mesh;
1608  const FiniteElement *fe;
1609  ElementTransformation *transf;
1610  FaceElementTransformations *face_elem_transf;
1611  Vector e_grad, a_grad, shape, el_dofs, err_val, ell_coeff_val;
1612  DenseMatrix dshape, dshapet, Jinv;
1613  Array<int> vdofs;
1614  IntegrationPoint eip;
1615  double error = 0.0;
1616 
1617  mesh = fes->GetMesh();
1618  dim = mesh->Dimension();
1619  e_grad.SetSize(dim);
1620  a_grad.SetSize(dim);
1621  Jinv.SetSize(dim);
1622 
1623  if (norm_type & 1)
1624  for (i = 0; i < mesh->GetNE(); i++)
1625  {
1626  fe = fes->GetFE(i);
1627  fdof = fe->GetDof();
1628  transf = mesh->GetElementTransformation(i);
1629  el_dofs.SetSize(fdof);
1630  dshape.SetSize(fdof, dim);
1631  dshapet.SetSize(fdof, dim);
1632  intorder = 2 * fe->GetOrder(); // <----------
1633  const IntegrationRule &ir = IntRules.Get(fe->GetGeomType(), intorder);
1634  fes->GetElementVDofs(i, vdofs);
1635  for (k = 0; k < fdof; k++)
1636  if (vdofs[k] >= 0)
1637  {
1638  el_dofs(k) = (*this)(vdofs[k]);
1639  }
1640  else
1641  {
1642  el_dofs(k) = - (*this)(-1-vdofs[k]);
1643  }
1644  for (j = 0; j < ir.GetNPoints(); j++)
1645  {
1646  const IntegrationPoint &ip = ir.IntPoint(j);
1647  fe->CalcDShape(ip, dshape);
1648  transf->SetIntPoint(&ip);
1649  exgrad->Eval(e_grad, *transf, ip);
1650  CalcInverse(transf->Jacobian(), Jinv);
1651  Mult(dshape, Jinv, dshapet);
1652  dshapet.MultTranspose(el_dofs, a_grad);
1653  e_grad -= a_grad;
1654  error += (ip.weight * transf->Weight() *
1655  ell_coeff->Eval(*transf, ip) *
1656  (e_grad * e_grad));
1657  }
1658  }
1659 
1660  if (norm_type & 2)
1661  for (i = 0; i < mesh->GetNFaces(); i++)
1662  {
1663  face_elem_transf = mesh->GetFaceElementTransformations(i, 5);
1664  int i1 = face_elem_transf->Elem1No;
1665  int i2 = face_elem_transf->Elem2No;
1666  intorder = fes->GetFE(i1)->GetOrder();
1667  if (i2 >= 0)
1668  if ( (k = fes->GetFE(i2)->GetOrder()) > intorder )
1669  {
1670  intorder = k;
1671  }
1672  intorder = 2 * intorder; // <-------------
1673  const IntegrationRule &ir =
1674  IntRules.Get(face_elem_transf->FaceGeom, intorder);
1675  err_val.SetSize(ir.GetNPoints());
1676  ell_coeff_val.SetSize(ir.GetNPoints());
1677  // side 1
1678  transf = face_elem_transf->Elem1;
1679  fe = fes->GetFE(i1);
1680  fdof = fe->GetDof();
1681  fes->GetElementVDofs(i1, vdofs);
1682  shape.SetSize(fdof);
1683  el_dofs.SetSize(fdof);
1684  for (k = 0; k < fdof; k++)
1685  if (vdofs[k] >= 0)
1686  {
1687  el_dofs(k) = (*this)(vdofs[k]);
1688  }
1689  else
1690  {
1691  el_dofs(k) = - (*this)(-1-vdofs[k]);
1692  }
1693  for (j = 0; j < ir.GetNPoints(); j++)
1694  {
1695  face_elem_transf->Loc1.Transform(ir.IntPoint(j), eip);
1696  fe->CalcShape(eip, shape);
1697  transf->SetIntPoint(&eip);
1698  ell_coeff_val(j) = ell_coeff->Eval(*transf, eip);
1699  err_val(j) = exsol->Eval(*transf, eip) - (shape * el_dofs);
1700  }
1701  if (i2 >= 0)
1702  {
1703  // side 2
1704  face_elem_transf = mesh->GetFaceElementTransformations(i, 10);
1705  transf = face_elem_transf->Elem2;
1706  fe = fes->GetFE(i2);
1707  fdof = fe->GetDof();
1708  fes->GetElementVDofs(i2, vdofs);
1709  shape.SetSize(fdof);
1710  el_dofs.SetSize(fdof);
1711  for (k = 0; k < fdof; k++)
1712  if (vdofs[k] >= 0)
1713  {
1714  el_dofs(k) = (*this)(vdofs[k]);
1715  }
1716  else
1717  {
1718  el_dofs(k) = - (*this)(-1-vdofs[k]);
1719  }
1720  for (j = 0; j < ir.GetNPoints(); j++)
1721  {
1722  face_elem_transf->Loc2.Transform(ir.IntPoint(j), eip);
1723  fe->CalcShape(eip, shape);
1724  transf->SetIntPoint(&eip);
1725  ell_coeff_val(j) += ell_coeff->Eval(*transf, eip);
1726  ell_coeff_val(j) *= 0.5;
1727  err_val(j) -= (exsol->Eval(*transf, eip) - (shape * el_dofs));
1728  }
1729  }
1730  face_elem_transf = mesh->GetFaceElementTransformations(i, 16);
1731  transf = face_elem_transf->Face;
1732  for (j = 0; j < ir.GetNPoints(); j++)
1733  {
1734  const IntegrationPoint &ip = ir.IntPoint(j);
1735  transf->SetIntPoint(&ip);
1736  error += (ip.weight * Nu * ell_coeff_val(j) *
1737  pow(transf->Weight(), 1.0-1.0/(dim-1)) *
1738  err_val(j) * err_val(j));
1739  }
1740  }
1741 
1742  if (error < 0.0)
1743  {
1744  return -sqrt(-error);
1745  }
1746  return sqrt(error);
1747 }
1748 
1750  Coefficient *exsol[], const IntegrationRule *irs[]) const
1751 {
1752  double error = 0.0, a;
1753  const FiniteElement *fe;
1754  ElementTransformation *transf;
1755  Vector shape;
1756  Array<int> vdofs;
1757  int fdof, d, i, intorder, j, k;
1758 
1759  for (i = 0; i < fes->GetNE(); i++)
1760  {
1761  fe = fes->GetFE(i);
1762  fdof = fe->GetDof();
1763  transf = fes->GetElementTransformation(i);
1764  shape.SetSize(fdof);
1765  intorder = 2*fe->GetOrder() + 1; // <----------
1766  const IntegrationRule *ir;
1767  if (irs)
1768  {
1769  ir = irs[fe->GetGeomType()];
1770  }
1771  else
1772  {
1773  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
1774  }
1775  fes->GetElementVDofs(i, vdofs);
1776  for (j = 0; j < ir->GetNPoints(); j++)
1777  {
1778  const IntegrationPoint &ip = ir->IntPoint(j);
1779  fe->CalcShape(ip, shape);
1780  transf->SetIntPoint(&ip);
1781  for (d = 0; d < fes->GetVDim(); d++)
1782  {
1783  a = 0;
1784  for (k = 0; k < fdof; k++)
1785  if (vdofs[fdof*d+k] >= 0)
1786  {
1787  a += (*this)(vdofs[fdof*d+k]) * shape(k);
1788  }
1789  else
1790  {
1791  a -= (*this)(-1-vdofs[fdof*d+k]) * shape(k);
1792  }
1793  a -= exsol[d]->Eval(*transf, ip);
1794  a = fabs(a);
1795  if (error < a)
1796  {
1797  error = a;
1798  }
1799  }
1800  }
1801  }
1802 
1803  return error;
1804 }
1805 
1807  Coefficient *exsol, VectorCoefficient *exgrad, int norm_type,
1808  Array<int> *elems, const IntegrationRule *irs[]) const
1809 {
1810  // assuming vdim is 1
1811  int i, fdof, dim, intorder, j, k;
1812  Mesh *mesh;
1813  const FiniteElement *fe;
1814  ElementTransformation *transf;
1815  Vector e_grad, a_grad, shape, el_dofs, err_val, ell_coeff_val;
1816  DenseMatrix dshape, dshapet, Jinv;
1817  Array<int> vdofs;
1818  double a, error = 0.0;
1819 
1820  mesh = fes->GetMesh();
1821  dim = mesh->Dimension();
1822  e_grad.SetSize(dim);
1823  a_grad.SetSize(dim);
1824  Jinv.SetSize(dim);
1825 
1826  if (norm_type & 1) // L_1 norm
1827  for (i = 0; i < mesh->GetNE(); i++)
1828  {
1829  if (elems != NULL && (*elems)[i] == 0) { continue; }
1830  fe = fes->GetFE(i);
1831  fdof = fe->GetDof();
1832  transf = fes->GetElementTransformation(i);
1833  el_dofs.SetSize(fdof);
1834  shape.SetSize(fdof);
1835  intorder = 2*fe->GetOrder() + 1; // <----------
1836  const IntegrationRule *ir;
1837  if (irs)
1838  {
1839  ir = irs[fe->GetGeomType()];
1840  }
1841  else
1842  {
1843  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
1844  }
1845  fes->GetElementVDofs(i, vdofs);
1846  for (k = 0; k < fdof; k++)
1847  if (vdofs[k] >= 0)
1848  {
1849  el_dofs(k) = (*this)(vdofs[k]);
1850  }
1851  else
1852  {
1853  el_dofs(k) = -(*this)(-1-vdofs[k]);
1854  }
1855  for (j = 0; j < ir->GetNPoints(); j++)
1856  {
1857  const IntegrationPoint &ip = ir->IntPoint(j);
1858  fe->CalcShape(ip, shape);
1859  transf->SetIntPoint(&ip);
1860  a = (el_dofs * shape) - (exsol->Eval(*transf, ip));
1861  error += ip.weight * transf->Weight() * fabs(a);
1862  }
1863  }
1864 
1865  if (norm_type & 2) // W^1_1 seminorm
1866  for (i = 0; i < mesh->GetNE(); i++)
1867  {
1868  if (elems != NULL && (*elems)[i] == 0) { continue; }
1869  fe = fes->GetFE(i);
1870  fdof = fe->GetDof();
1871  transf = mesh->GetElementTransformation(i);
1872  el_dofs.SetSize(fdof);
1873  dshape.SetSize(fdof, dim);
1874  dshapet.SetSize(fdof, dim);
1875  intorder = 2*fe->GetOrder() + 1; // <----------
1876  const IntegrationRule *ir;
1877  if (irs)
1878  {
1879  ir = irs[fe->GetGeomType()];
1880  }
1881  else
1882  {
1883  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
1884  }
1885  fes->GetElementVDofs(i, vdofs);
1886  for (k = 0; k < fdof; k++)
1887  if (vdofs[k] >= 0)
1888  {
1889  el_dofs(k) = (*this)(vdofs[k]);
1890  }
1891  else
1892  {
1893  el_dofs(k) = -(*this)(-1-vdofs[k]);
1894  }
1895  for (j = 0; j < ir->GetNPoints(); j++)
1896  {
1897  const IntegrationPoint &ip = ir->IntPoint(j);
1898  fe->CalcDShape(ip, dshape);
1899  transf->SetIntPoint(&ip);
1900  exgrad->Eval(e_grad, *transf, ip);
1901  CalcInverse(transf->Jacobian(), Jinv);
1902  Mult(dshape, Jinv, dshapet);
1903  dshapet.MultTranspose(el_dofs, a_grad);
1904  e_grad -= a_grad;
1905  error += ip.weight * transf->Weight() * e_grad.Norml1();
1906  }
1907  }
1908 
1909  return error;
1910 }
1911 
1912 double GridFunction::ComputeLpError(const double p, Coefficient &exsol,
1913  Coefficient *weight,
1914  const IntegrationRule *irs[]) const
1915 {
1916  double error = 0.0;
1917  const FiniteElement *fe;
1919  Vector vals;
1920 
1921  for (int i = 0; i < fes->GetNE(); i++)
1922  {
1923  fe = fes->GetFE(i);
1924  const IntegrationRule *ir;
1925  if (irs)
1926  {
1927  ir = irs[fe->GetGeomType()];
1928  }
1929  else
1930  {
1931  int intorder = 2*fe->GetOrder() + 1; // <----------
1932  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
1933  }
1934  GetValues(i, *ir, vals);
1935  T = fes->GetElementTransformation(i);
1936  for (int j = 0; j < ir->GetNPoints(); j++)
1937  {
1938  const IntegrationPoint &ip = ir->IntPoint(j);
1939  T->SetIntPoint(&ip);
1940  double err = fabs(vals(j) - exsol.Eval(*T, ip));
1941  if (p < numeric_limits<double>::infinity())
1942  {
1943  err = pow(err, p);
1944  if (weight)
1945  {
1946  err *= weight->Eval(*T, ip);
1947  }
1948  error += ip.weight * T->Weight() * err;
1949  }
1950  else
1951  {
1952  if (weight)
1953  {
1954  err *= weight->Eval(*T, ip);
1955  }
1956  error = std::max(error, err);
1957  }
1958  }
1959  }
1960 
1961  if (p < numeric_limits<double>::infinity())
1962  {
1963  // negative quadrature weights may cause the error to be negative
1964  if (error < 0.)
1965  {
1966  error = -pow(-error, 1./p);
1967  }
1968  else
1969  {
1970  error = pow(error, 1./p);
1971  }
1972  }
1973 
1974  return error;
1975 }
1976 
1977 double GridFunction::ComputeLpError(const double p, VectorCoefficient &exsol,
1978  Coefficient *weight,
1979  VectorCoefficient *v_weight,
1980  const IntegrationRule *irs[]) const
1981 {
1982  double error = 0.0;
1983  const FiniteElement *fe;
1985  DenseMatrix vals, exact_vals;
1986  Vector loc_errs;
1987 
1988  for (int i = 0; i < fes->GetNE(); i++)
1989  {
1990  fe = fes->GetFE(i);
1991  const IntegrationRule *ir;
1992  if (irs)
1993  {
1994  ir = irs[fe->GetGeomType()];
1995  }
1996  else
1997  {
1998  int intorder = 2*fe->GetOrder() + 1; // <----------
1999  ir = &(IntRules.Get(fe->GetGeomType(), intorder));
2000  }
2001  T = fes->GetElementTransformation(i);
2002  GetVectorValues(*T, *ir, vals);
2003  exsol.Eval(exact_vals, *T, *ir);
2004  vals -= exact_vals;
2005  loc_errs.SetSize(vals.Width());
2006  if (!v_weight)
2007  {
2008  // compute the lengths of the errors at the integration points
2009  // thus the vector norm is rotationally invariant
2010  vals.Norm2(loc_errs);
2011  }
2012  else
2013  {
2014  v_weight->Eval(exact_vals, *T, *ir);
2015  // column-wise dot product of the vector error (in vals) and the
2016  // vector weight (in exact_vals)
2017  for (int j = 0; j < vals.Width(); j++)
2018  {
2019  double err = 0.0;
2020  for (int d = 0; d < vals.Height(); d++)
2021  {
2022  err += vals(d,j)*exact_vals(d,j);
2023  }
2024  loc_errs(j) = fabs(err);
2025  }
2026  }
2027  for (int j = 0; j < ir->GetNPoints(); j++)
2028  {
2029  const IntegrationPoint &ip = ir->IntPoint(j);
2030  T->SetIntPoint(&ip);
2031  double err = loc_errs(j);
2032  if (p < numeric_limits<double>::infinity())
2033  {
2034  err = pow(err, p);
2035  if (weight)
2036  {
2037  err *= weight->Eval(*T, ip);
2038  }
2039  error += ip.weight * T->Weight() * err;
2040  }
2041  else
2042  {
2043  if (weight)
2044  {
2045  err *= weight->Eval(*T, ip);
2046  }
2047  error = std::max(error, err);
2048  }
2049  }
2050  }
2051 
2052  if (p < numeric_limits<double>::infinity())
2053  {
2054  // negative quadrature weights may cause the error to be negative
2055  if (error < 0.)
2056  {
2057  error = -pow(-error, 1./p);
2058  }
2059  else
2060  {
2061  error = pow(error, 1./p);
2062  }
2063  }
2064 
2065  return error;
2066 }
2067 
2069 {
2070  for (int i = 0; i < size; i++)
2071  {
2072  data[i] = value;
2073  }
2074  return *this;
2075 }
2076 
2078 {
2079  for (int i = 0; i < size; i++)
2080  {
2081  data[i] = v(i);
2082  }
2083  return *this;
2084 }
2085 
2087 {
2088  return this->operator=((const Vector &)v);
2089 }
2090 
2092 {
2094  if (P)
2095  {
2096  this->SetSize(P->Height());
2097  P->Mult(x, *this);
2098  }
2099  else // assume conforming mesh
2100  {
2101  *this = x;
2102  }
2103 }
2104 
2106 {
2107  if (fes->Nonconforming())
2108  {
2109  Vector x = *this;
2111  }
2112 }
2113 
2115 {
2117  if (R)
2118  {
2119  x.SetSize(R->Height());
2120  R->Mult(*this, x);
2121  }
2122  else // assume conforming mesh
2123  {
2124  x = *this;
2125  }
2126 }
2127 
2129 {
2130  if (fes->Nonconforming())
2131  {
2132  Vector x;
2133  ConformingProject(x);
2134  static_cast<Vector&>(*this) = x;
2135  }
2136 }
2137 
2138 void GridFunction::Save(std::ostream &out) const
2139 {
2140  fes->Save(out);
2141  out << '\n';
2142  if (fes->GetOrdering() == Ordering::byNODES)
2143  {
2144  Vector::Print(out, 1);
2145  }
2146  else
2147  {
2148  Vector::Print(out, fes->GetVDim());
2149  }
2150 }
2151 
2152 void GridFunction::SaveVTK(std::ostream &out, const std::string &field_name,
2153  int ref)
2154 {
2155  Mesh *mesh = fes->GetMesh();
2156  RefinedGeometry *RefG;
2157  Vector val;
2158  DenseMatrix vval, pmat;
2159  int vec_dim = VectorDim();
2160 
2161  if (vec_dim == 1)
2162  {
2163  // scalar data
2164  out << "SCALARS " << field_name << " double 1\n"
2165  << "LOOKUP_TABLE default\n";
2166  for (int i = 0; i < mesh->GetNE(); i++)
2167  {
2168  RefG = GlobGeometryRefiner.Refine(
2169  mesh->GetElementBaseGeometry(i), ref, 1);
2170 
2171  GetValues(i, RefG->RefPts, val, pmat);
2172 
2173  for (int j = 0; j < val.Size(); j++)
2174  {
2175  out << val(j) << '\n';
2176  }
2177  }
2178  }
2179  else if (vec_dim == mesh->Dimension())
2180  {
2181  // vector data
2182  out << "VECTORS " << field_name << " double\n";
2183  for (int i = 0; i < mesh->GetNE(); i++)
2184  {
2185  RefG = GlobGeometryRefiner.Refine(
2186  mesh->GetElementBaseGeometry(i), ref, 1);
2187 
2188  GetVectorValues(i, RefG->RefPts, vval, pmat);
2189 
2190  for (int j = 0; j < vval.Width(); j++)
2191  {
2192  out << vval(0, j) << ' ' << vval(1, j) << ' ';
2193  if (vval.Height() == 2)
2194  {
2195  out << 0.0;
2196  }
2197  else
2198  {
2199  out << vval(2, j);
2200  }
2201  out << '\n';
2202  }
2203  }
2204  }
2205  else
2206  {
2207  // other data: save the components as separate scalars
2208  for (int vd = 0; vd < vec_dim; vd++)
2209  {
2210  out << "SCALARS " << field_name << vd << " double 1\n"
2211  << "LOOKUP_TABLE default\n";
2212  for (int i = 0; i < mesh->GetNE(); i++)
2213  {
2214  RefG = GlobGeometryRefiner.Refine(
2215  mesh->GetElementBaseGeometry(i), ref, 1);
2216 
2217  GetValues(i, RefG->RefPts, val, pmat, vd + 1);
2218 
2219  for (int j = 0; j < val.Size(); j++)
2220  {
2221  out << val(j) << '\n';
2222  }
2223  }
2224  }
2225  }
2226 }
2227 
2228 void GridFunction::SaveSTLTri(std::ostream &out, double p1[], double p2[],
2229  double p3[])
2230 {
2231  double v1[3] = { p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2] };
2232  double v2[3] = { p3[0] - p1[0], p3[1] - p1[1], p3[2] - p1[2] };
2233  double n[] = { v1[1] * v2[2] - v1[2] * v2[1],
2234  v1[2] * v2[0] - v1[0] * v2[2],
2235  v1[0] * v2[1] - v1[1] * v2[0]
2236  };
2237  double rl = 1.0 / sqrt(n[0] * n[0] + n[1] * n[1] + n[2] * n[2]);
2238  n[0] *= rl; n[1] *= rl; n[2] *= rl;
2239 
2240  out << " facet normal " << n[0] << ' ' << n[1] << ' ' << n[2]
2241  << "\n outer loop"
2242  << "\n vertex " << p1[0] << ' ' << p1[1] << ' ' << p1[2]
2243  << "\n vertex " << p2[0] << ' ' << p2[1] << ' ' << p2[2]
2244  << "\n vertex " << p3[0] << ' ' << p3[1] << ' ' << p3[2]
2245  << "\n endloop\n endfacet\n";
2246 }
2247 
2248 void GridFunction::SaveSTL(std::ostream &out, int TimesToRefine)
2249 {
2250  Mesh *mesh = fes->GetMesh();
2251 
2252  if (mesh->Dimension() != 2)
2253  {
2254  return;
2255  }
2256 
2257  int i, j, k, l, n;
2258  DenseMatrix pointmat;
2259  Vector values;
2260  RefinedGeometry * RefG;
2261  double pts[4][3], bbox[3][2];
2262 
2263  out << "solid GridFunction\n";
2264 
2265  bbox[0][0] = bbox[0][1] = bbox[1][0] = bbox[1][1] =
2266  bbox[2][0] = bbox[2][1] = 0.0;
2267  for (i = 0; i < mesh->GetNE(); i++)
2268  {
2269  n = fes->GetFE(i)->GetGeomType();
2270  RefG = GlobGeometryRefiner.Refine(n, TimesToRefine);
2271  GetValues(i, RefG->RefPts, values, pointmat);
2272  Array<int> &RG = RefG->RefGeoms;
2273  n = Geometries.NumBdr(n);
2274  for (k = 0; k < RG.Size()/n; k++)
2275  {
2276  for (j = 0; j < n; j++)
2277  {
2278  l = RG[n*k+j];
2279  pts[j][0] = pointmat(0,l);
2280  pts[j][1] = pointmat(1,l);
2281  pts[j][2] = values(l);
2282  }
2283 
2284  if (n == 3)
2285  {
2286  SaveSTLTri(out, pts[0], pts[1], pts[2]);
2287  }
2288  else
2289  {
2290  SaveSTLTri(out, pts[0], pts[1], pts[2]);
2291  SaveSTLTri(out, pts[0], pts[2], pts[3]);
2292  }
2293  }
2294 
2295  if (i == 0)
2296  {
2297  bbox[0][0] = pointmat(0,0);
2298  bbox[0][1] = pointmat(0,0);
2299  bbox[1][0] = pointmat(1,0);
2300  bbox[1][1] = pointmat(1,0);
2301  bbox[2][0] = values(0);
2302  bbox[2][1] = values(0);
2303  }
2304 
2305  for (j = 0; j < values.Size(); j++)
2306  {
2307  if (bbox[0][0] > pointmat(0,j))
2308  {
2309  bbox[0][0] = pointmat(0,j);
2310  }
2311  if (bbox[0][1] < pointmat(0,j))
2312  {
2313  bbox[0][1] = pointmat(0,j);
2314  }
2315  if (bbox[1][0] > pointmat(1,j))
2316  {
2317  bbox[1][0] = pointmat(1,j);
2318  }
2319  if (bbox[1][1] < pointmat(1,j))
2320  {
2321  bbox[1][1] = pointmat(1,j);
2322  }
2323  if (bbox[2][0] > values(j))
2324  {
2325  bbox[2][0] = values(j);
2326  }
2327  if (bbox[2][1] < values(j))
2328  {
2329  bbox[2][1] = values(j);
2330  }
2331  }
2332  }
2333 
2334  cout << "[xmin,xmax] = [" << bbox[0][0] << ',' << bbox[0][1] << "]\n"
2335  << "[ymin,ymax] = [" << bbox[1][0] << ',' << bbox[1][1] << "]\n"
2336  << "[zmin,zmax] = [" << bbox[2][0] << ',' << bbox[2][1] << ']'
2337  << endl;
2338 
2339  out << "endsolid GridFunction" << endl;
2340 }
2341 
2342 std::ostream &operator<<(std::ostream &out, const GridFunction &sol)
2343 {
2344  sol.Save(out);
2345  return out;
2346 }
2347 
2348 
2349 
2351  GridFunction &u,
2352  GridFunction &flux, Vector &error_estimates,
2353  Array<int>* aniso_flags,
2354  int with_subdomains)
2355 {
2356  FiniteElementSpace *ufes = u.FESpace();
2357  FiniteElementSpace *ffes = flux.FESpace();
2358  ElementTransformation *Transf;
2359 
2360  int dim = ufes->GetMesh()->Dimension();
2361  int nfe = ufes->GetNE();
2362 
2363  Array<int> udofs;
2364  Array<int> fdofs;
2365  Vector ul, fl, fla, d_xyz;
2366 
2367  error_estimates.SetSize(nfe);
2368  if (aniso_flags)
2369  {
2370  aniso_flags->SetSize(nfe);
2371  d_xyz.SetSize(dim);
2372  }
2373 
2374  int nsd = 1;
2375  if (with_subdomains)
2376  {
2377  for (int i = 0; i < nfe; i++)
2378  {
2379  int attr = ufes->GetAttribute(i);
2380  if (attr > nsd) { nsd = attr; }
2381  }
2382  }
2383 
2384  for (int s = 1; s <= nsd; s++)
2385  {
2386  // This calls the parallel version when u is a ParGridFunction
2387  u.ComputeFlux(blfi, flux, 0, (with_subdomains ? s : 0));
2388 
2389  for (int i = 0; i < nfe; i++)
2390  {
2391  if (with_subdomains && ufes->GetAttribute(i) != s) { continue; }
2392 
2393  ufes->GetElementVDofs(i, udofs);
2394  ffes->GetElementVDofs(i, fdofs);
2395 
2396  u.GetSubVector(udofs, ul);
2397  flux.GetSubVector(fdofs, fla);
2398 
2399  Transf = ufes->GetElementTransformation(i);
2400  blfi.ComputeElementFlux(*ufes->GetFE(i), *Transf, ul,
2401  *ffes->GetFE(i), fl, 0);
2402 
2403  fl -= fla;
2404 
2405  error_estimates(i) =
2406  blfi.ComputeFluxEnergy(*ffes->GetFE(i), *Transf, fl,
2407  (aniso_flags ? &d_xyz : NULL));
2408 
2409  if (aniso_flags)
2410  {
2411  double sum = 0;
2412  for (int k = 0; k < dim; k++)
2413  {
2414  sum += d_xyz[k];
2415  }
2416 
2417  double thresh = 0.15 * 3.0/dim;
2418  int flag = 0;
2419  for (int k = 0; k < dim; k++)
2420  {
2421  if (d_xyz[k] / sum > thresh) { flag |= (1 << k); }
2422  }
2423 
2424  (*aniso_flags)[i] = flag;
2425  }
2426  }
2427  }
2428 }
2429 
2430 
2431 double ComputeElementLpDistance(double p, int i,
2432  GridFunction& gf1, GridFunction& gf2)
2433 {
2434  double norm = 0.0;
2435 
2436  FiniteElementSpace *fes1 = gf1.FESpace();
2437  FiniteElementSpace *fes2 = gf2.FESpace();
2438 
2439  const FiniteElement* fe1 = fes1->GetFE(i);
2440  const FiniteElement* fe2 = fes2->GetFE(i);
2441 
2442  const IntegrationRule *ir;
2443  int intorder = 2*std::max(fe1->GetOrder(),fe2->GetOrder()) + 1; // <-------
2444  ir = &(IntRules.Get(fe1->GetGeomType(), intorder));
2445  int nip = ir->GetNPoints();
2446  Vector val1, val2;
2447 
2448 
2450  for (int j = 0; j < nip; j++)
2451  {
2452  const IntegrationPoint &ip = ir->IntPoint(j);
2453  T->SetIntPoint(&ip);
2454 
2455  gf1.GetVectorValue(i, ip, val1);
2456  gf2.GetVectorValue(i, ip, val2);
2457 
2458  val1 -= val2;
2459  double err = val1.Norml2();
2460  if (p < numeric_limits<double>::infinity())
2461  {
2462  err = pow(err, p);
2463  norm += ip.weight * T->Weight() * err;
2464  }
2465  else
2466  {
2467  norm = std::max(norm, err);
2468  }
2469  }
2470 
2471  if (p < numeric_limits<double>::infinity())
2472  {
2473  // Negative quadrature weights may cause the norm to be negative
2474  if (norm < 0.)
2475  {
2476  norm = -pow(-norm, 1./p);
2477  }
2478  else
2479  {
2480  norm = pow(norm, 1./p);
2481  }
2482  }
2483 
2484  return norm;
2485 }
2486 
2487 
2489  const IntegrationPoint &ip)
2490 {
2491  ElementTransformation *T_in =
2492  mesh_in->GetElementTransformation(T.ElementNo / n);
2493  T_in->SetIntPoint(&ip);
2494  return sol_in.Eval(*T_in, ip);
2495 }
2496 
2497 
2499  GridFunction *sol, const int ny)
2500 {
2501  GridFunction *sol2d;
2502 
2503  FiniteElementCollection *solfec2d;
2504  const char *name = sol->FESpace()->FEColl()->Name();
2505  string cname = name;
2506  if (cname == "Linear")
2507  {
2508  solfec2d = new LinearFECollection;
2509  }
2510  else if (cname == "Quadratic")
2511  {
2512  solfec2d = new QuadraticFECollection;
2513  }
2514  else if (cname == "Cubic")
2515  {
2516  solfec2d = new CubicFECollection;
2517  }
2518  else if (!strncmp(name, "H1_", 3))
2519  {
2520  solfec2d = new H1_FECollection(atoi(name + 7), 2);
2521  }
2522  else if (!strncmp(name, "L2_T", 4))
2523  {
2524  solfec2d = new L2_FECollection(atoi(name + 10), 2);
2525  }
2526  else if (!strncmp(name, "L2_", 3))
2527  {
2528  solfec2d = new L2_FECollection(atoi(name + 7), 2);
2529  }
2530  else
2531  {
2532  cerr << "Extrude1DGridFunction : unknown FE collection : "
2533  << cname << endl;
2534  return NULL;
2535  }
2536  FiniteElementSpace *solfes2d;
2537  // assuming sol is scalar
2538  solfes2d = new FiniteElementSpace(mesh2d, solfec2d);
2539  sol2d = new GridFunction(solfes2d);
2540  sol2d->MakeOwner(solfec2d);
2541  {
2542  GridFunctionCoefficient csol(sol);
2543  ExtrudeCoefficient c2d(mesh, csol, ny);
2544  sol2d->ProjectCoefficient(c2d);
2545  }
2546  return sol2d;
2547 }
2548 
2549 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:228
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
Definition: gridfunc.cpp:2228
Abstract class for Finite Elements.
Definition: fe.hpp:44
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:9824
int Size() const
Logical size of the array.
Definition: array.hpp:109
int GetVSize() const
Definition: fespace.hpp:164
int GetDim() const
Returns the space dimension for the finite element.
Definition: fe.hpp:82
int GetAttribute(int i) const
Definition: fespace.hpp:218
ElementTransformation * Face
Definition: eltrans.hpp:119
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:162
Class for integration rule.
Definition: intrules.hpp:83
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:479
void NewDataAndSize(double *d, int s)
Definition: vector.hpp:72
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:116
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:30
Class used for extruding scalar GridFunctions.
Definition: gridfunc.hpp:312
void GetElementVertices(int i, Array< int > &vertices) const
Returns the vertices of element i.
Definition: fespace.hpp:198
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:235
void MergeGridFunctions(GridFunction *gf_array[], int num_pieces, GridFunction &merged)
Definition: nurbs.cpp:1568
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:1077
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:75
void Print(std::ostream &out=std::cout, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:583
void ConformingProlongate()
Definition: gridfunc.cpp:2105
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:292
void SetSize(int s)
Resizes the vector if the new size is different.
Definition: vector.hpp:259
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:161
GridFunction * Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d, GridFunction *sol, const int ny)
Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
Definition: gridfunc.cpp:2498
void MakeOwner(FiniteElementCollection *_fec)
Make the GridFunction the owner of &#39;fec&#39; and &#39;fes&#39;.
Definition: gridfunc.hpp:70
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh)
Definition: gridfunc.cpp:774
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: gridfunc.cpp:209
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1269
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:445
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:644
int GetElementBaseGeometry(int i) const
Definition: mesh.hpp:494
void GetBdrValuesFrom(GridFunction &)
Definition: gridfunc.cpp:537
void GetGradient(ElementTransformation &tr, Vector &grad)
Definition: gridfunc.cpp:876
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:33
RefinedGeometry * Refine(int Geom, int Times, int ETimes=1)
Definition: geom.cpp:465
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:457
void StealData(double **p)
Changes the ownership of the data; after the call the Vector is empty.
Definition: vector.hpp:97
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:91
int VectorDim() const
Definition: gridfunc.cpp:224
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:260
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
Delta function coefficient.
void ImposeBounds(int i, const Vector &weights, const Vector &_lo, const Vector &_hi)
Definition: gridfunc.cpp:993
int Size() const
Returns the size of the vector.
Definition: vector.hpp:84
int GetBdrAttribute(int i) const
Definition: fespace.hpp:220
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:120
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:454
void GetElementVertices(int i, Array< int > &dofs) const
Returns the indices of the dofs of element i.
Definition: mesh.hpp:501
int GetNV() const
Returns number of nodes in the mesh.
Definition: fespace.hpp:188
const IntegrationPoint & GetIntPoint()
Definition: eltrans.hpp:35
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Definition: bilininteg.hpp:74
ElementTransformation * Elem2
Definition: eltrans.hpp:119
double * GetData() const
Definition: vector.hpp:88
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:3100
void GetNodalValues(int i, Array< double > &nval, int vdim=1) const
Returns the values in the vertices of i&#39;th element for dimension vdim.
Definition: gridfunc.cpp:240
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:557
const IntegrationRule * GetVertices(int GeomType)
Definition: geom.cpp:172
void GetEdgeTransformation(int i, IsoparametricTransformation *EdTr)
Definition: mesh.cpp:369
void GetVectorFieldNodalValues(Vector &val, int comp) const
Definition: gridfunc.cpp:643
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2138
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Definition: vector.cpp:45
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:327
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:185
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Definition: gridfunc.cpp:2488
double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:205
Geometry Geometries
Definition: geom.cpp:443
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:231
void GetDerivative(int comp, int der_comp, GridFunction &der)
Definition: gridfunc.cpp:714
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:193
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1383
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
Definition: gridfunc.cpp:162
void SetLinearConstraint(const Vector &_w, double _a)
Definition: solvers.cpp:1413
int dim
Definition: ex3.cpp:48
double GetDivergence(ElementTransformation &tr)
Definition: gridfunc.cpp:793
const FiniteElement * GetEdgeElement(int i) const
Definition: fespace.cpp:1468
virtual ~GridFunction()
Destroys grid function.
Definition: gridfunc.cpp:128
IntegrationRules IntRules(0)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:292
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
const SparseMatrix * GetConformingRestriction()
Definition: fespace.cpp:919
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:190
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:191
Data type sparse matrix.
Definition: sparsemat.hpp:38
const double * Center()
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:140
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
Definition: gridfunc.cpp:574
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:120
Coefficient * Weight()
double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:1601
void SetMaxIter(int max_it)
Definition: solvers.hpp:61
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
void Norm2(double *v) const
Take the 2-norm of the columns of A and store in v.
Definition: densemat.cpp:678
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:1100
const IntegrationRule & GetNodes() const
Definition: fe.hpp:113
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
Definition: bilininteg.hpp:68
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:914
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i&#39;th boundary element.
Definition: fespace.hpp:215
IntegrationRule RefPts
Definition: geom.hpp:80
void GetCurl(ElementTransformation &tr, Vector &curl)
Definition: gridfunc.cpp:825
const SparseMatrix * GetConformingProlongation()
Definition: fespace.cpp:912
virtual double Weight()=0
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:324
int GetVDim()
Returns dimension of the vector.
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:154
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:63
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:270
int Dimension() const
Definition: mesh.hpp:475
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
Definition: ncmesh.cpp:2387
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:527
double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:1912
void GetEdgeDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1354
void ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains)
Definition: gridfunc.cpp:2350
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2972
int GetGeomType() const
Returns the geometry type:
Definition: fe.hpp:85
int GetLocalDofForDof(int i)
Definition: fespace.hpp:282
int GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:176
int SpaceDimension() const
Definition: mesh.hpp:476
int GetElementForDof(int i)
Definition: fespace.hpp:281
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:713
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i&#39;th element.
Definition: fespace.hpp:206
double Norml1() const
Returns the l_1 norm of the vector.
Definition: vector.cpp:662
bool FaceIsInterior(int FaceNo) const
Return true if the given face is interior.
Definition: mesh.hpp:633
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
void SetAbsTol(double atol)
Definition: solvers.hpp:60
void SetRelTol(double rtol)
Definition: solvers.hpp:59
double Distance(const double *x, const double *y, const int n)
Definition: vector.hpp:325
Abstract finite element space.
Definition: fespace.hpp:62
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:417
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
FiniteElementSpace * fes
FE space on which grid function lives.
Definition: gridfunc.hpp:31
void ProjectGridFunction(const GridFunction &src)
Definition: gridfunc.cpp:960
void mfem_error(const char *msg)
Definition: error.cpp:23
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:323
void SetSubVector(const Array< int > &dofs, const Vector &elemvect)
Definition: vector.cpp:493
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:243
int NumBdr(int GeomType)
Definition: geom.hpp:71
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:200
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:44
FiniteElementCollection * fec
Used when the grid function is read from a file.
Definition: gridfunc.hpp:34
double ComputeElementLpDistance(double p, int i, GridFunction &gf1, GridFunction &gf2)
Compute the Lp distance between two grid functions on the given element.
Definition: gridfunc.cpp:2431
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:457
static FiniteElementCollection * New(const char *name)
Definition: fe_coll.cpp:44
void SetBounds(const Vector &_lo, const Vector &_hi)
Definition: solvers.cpp:1407
bool Nonconforming() const
Definition: fespace.hpp:146
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1212
NURBSExtension * NURBSext
Definition: mesh.hpp:142
int GetNV() const
Definition: mesh.hpp:451
int GetNVDofs() const
Definition: fespace.hpp:180
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:51
virtual const char * Name() const
Definition: fe_coll.hpp:36
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:200
Class for integration point with weight.
Definition: intrules.hpp:25
void SaveSTL(std::ostream &out, int TimesToRefine=1)
Definition: gridfunc.cpp:2248
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:173
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:700
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:531
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition: gridfunc.cpp:355
void ConformingProject()
Definition: gridfunc.cpp:2128
void GetElementAverages(GridFunction &avgs)
Definition: gridfunc.cpp:930
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1199
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1454
ElementTransformation * Elem1
Definition: eltrans.hpp:119
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
Definition: gridfunc.cpp:405
virtual void Mult(const Vector &xt, Vector &x) const
Operator application.
Definition: solvers.cpp:1438
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1140
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:283
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Definition: gridfunc.cpp:616
NCMesh * ncmesh
Definition: mesh.hpp:143
double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:1806
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:179
GridFunction & operator=(double value)
Redefine &#39;=&#39; for GridFunction = constant.
Definition: gridfunc.cpp:2068
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:463
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:3362
void DofsToVDofs(Array< int > &dofs) const
Definition: fespace.cpp:35
Vector data type.
Definition: vector.hpp:33
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:178
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:142
virtual void Transform(const IntegrationPoint &, Vector &)=0
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:54
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad)
Definition: gridfunc.cpp:894
virtual const DenseMatrix & Jacobian()=0
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1419
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:167
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:71
double * data
Definition: vector.hpp:38
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Definition: gridfunc.cpp:918
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
int GetRangeType() const
Definition: fe.hpp:96
void Save(std::ostream &out) const
Definition: fespace.cpp:1635
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:726
void GetValuesFrom(GridFunction &)
Definition: gridfunc.cpp:500
int GetNFDofs() const
Definition: fespace.hpp:182
int GetNEDofs() const
Definition: fespace.hpp:181
Array< int > RefGeoms
Definition: geom.hpp:81
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:148
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:87
void SaveVTK(std::ostream &out, const std::string &field_name, int ref)
Definition: gridfunc.cpp:2152
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Definition: gridfunc.hpp:170
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:279
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:95
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
Definition: gridfunc.cpp:674