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