MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
display-basis.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 // ----------------------------------------------------------------
13 // Display Basis Miniapp: Visualize finite element basis functions
14 // ----------------------------------------------------------------
15 //
16 // This miniapp visualizes various types of finite element basis functions on a
17 // single mesh element in 1D, 2D and 3D. The order and the type of finite
18 // element space can be changed, and the mesh element is either the reference
19 // one, or a simple transformation of it. Dynamic creation and interaction with
20 // multiple GLVis windows is demonstrated.
21 //
22 // Compile with: make display-basis
23 //
24 // Sample runs: display-basis
25 // display_basis -e 2 -b 3 -o 3
26 // display-basis -e 5 -b 1 -o 1
27 
28 #include "mfem.hpp"
29 #include "../common/fem_extras.hpp"
30 #include "../common/mesh_extras.hpp"
31 #include <vector>
32 #include <iostream>
33 
34 using namespace std;
35 using namespace mfem;
36 using namespace mfem::miniapps;
37 
38 // Data structure used to collect visualization window layout parameters
39 struct VisWinLayout
40 {
41  int nx;
42  int ny;
43  int w;
44  int h;
45 };
46 
47 // Data structure used to define simple coordinate transformations
48 struct DeformationData
49 {
50  double uniformScale;
51 
52  int squeezeAxis;
53  double squeezeFactor;
54 
55  int shearAxis;
56  Vector shearVec;
57 };
58 
59 /** The Deformation class implements three simple coordinate transformations:
60  Uniform Scaling:
61  u = a v for a scalar constant 'a'
62 
63  Compression or Squeeze (along a coordinate axis):
64  / 1/b 0 \ / 1/b 0 0 \ for a scalar constant b
65  u = \ 0 b / v or u = | 0 c 0 | v, and c = sqrt(b)
66  \ 0 0 c / the axis can also be chosen
67 
68  Shear:
69  u = v + v_i * s where 's' is the shear vector
70  and 'i' is the shear axis
71 */
72 class Deformation : public VectorCoefficient
73 {
74 public:
75 
76  enum DefType {INVALID, UNIFORM, SQUEEZE, SHEAR};
77 
78  Deformation(int dim, DefType dType, const DeformationData & data)
79  : VectorCoefficient(dim), dim_(dim), dType_(dType), data_(data) {}
80 
81  void Eval(Vector &v, ElementTransformation &T, const IntegrationPoint &ip);
82  using VectorCoefficient::Eval;
83 private:
84  void Def1D(const Vector & u, Vector & v);
85  void Def2D(const Vector & u, Vector & v);
86  void Def3D(const Vector & u, Vector & v);
87 
88  int dim_;
89  DefType dType_;
90  const DeformationData & data_;
91 };
92 
93 string elemTypeStr(const Element::Type & eType);
94 inline bool elemIs1D(const Element::Type & eType);
95 inline bool elemIs2D(const Element::Type & eType);
96 inline bool elemIs3D(const Element::Type & eType);
97 
98 string basisTypeStr(char bType);
99 inline bool basisIs1D(char bType);
100 inline bool basisIs2D(char bType);
101 inline bool basisIs3D(char bType);
102 
103 string mapTypeStr(int mType);
104 
105 int update_basis(vector<socketstream*> & sock, const VisWinLayout & vwl,
106  Element::Type e, char bType, int bOrder, int mType,
107  Deformation::DefType dType, const DeformationData & defData,
108  bool visualization);
109 
110 int main(int argc, char *argv[])
111 {
112  // Parse command-line options.
113  Element::Type eType = Element::TRIANGLE;
114  char bType = 'h';
115  int bOrder = 2;
116  int mType = 0;
117 
118  int eInt = -1;
119  int bInt = -1;
120 
121  VisWinLayout vwl;
122  vwl.nx = 5;
123  vwl.ny = 3;
124  vwl.w = 250;
125  vwl.h = 250;
126 
127  Deformation::DefType dType = Deformation::INVALID;
128  DeformationData defData;
129 
130  bool visualization = true;
131 
132  vector<socketstream*> sock;
133 
134  OptionsParser args(argc, argv);
135  args.AddOption(&eInt, "-e", "--elem-type",
136  "Element Type: (1-Segment, 2-Triangle, 3-Quadrilateral, "
137  "4-Tetrahedron, 5-Hexahedron)");
138  args.AddOption(&bInt, "-b", "--basis-type",
139  "Basis Function Type (0-H1, 1-Nedelec, 2-Raviart-Thomas, "
140  "3-L2, 4-Fixed Order Cont.,\n\t5-Gaussian Discontinuous (2D),"
141  " 6-Crouzeix-Raviart)");
142  args.AddOption(&bOrder, "-o", "--order", "Basis function order");
143  args.AddOption(&vwl.nx, "-nx", "--num-win-x",
144  "Number of Viz windows in X");
145  args.AddOption(&vwl.ny, "-ny", "--num-win-y",
146  "Number of Viz windows in y");
147  args.AddOption(&vwl.w, "-w", "--width",
148  "Width of Viz windows");
149  args.AddOption(&vwl.h, "-h", "--height",
150  "Height of Viz windows");
151  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
152  "--no-visualization",
153  "Enable or disable GLVis visualization.");
154  args.Parse();
155  if (!args.Good())
156  {
157  args.PrintUsage(cout);
158  return 1;
159  }
160  {
161  args.PrintOptions(cout);
162  }
163  if ( eInt > 0 && eInt < 6 )
164  {
165  eType = (Element::Type)eInt;
166  }
167  switch (bInt)
168  {
169  case 0:
170  bType = 'h';
171  break;
172  case 1:
173  bType = 'n';
174  break;
175  case 2:
176  bType = 'r';
177  break;
178  case 3:
179  bType = 'l';
180  break;
181  case 4:
182  bType = 'f';
183  break;
184  case 5:
185  bType = 'g';
186  break;
187  case 6:
188  bType = 'c';
189  break;
190  default:
191  bType = 'h';
192  }
193 
194  // Collect user input
195  bool print_char = true;
196  while (true)
197  {
198  if (print_char)
199  {
200  cout << endl;
201  cout << "Element Type: " << elemTypeStr(eType) << endl;
202  cout << "Basis Type: " << basisTypeStr(bType) << endl;;
203  cout << "Basis function order: " << bOrder << endl;
204  cout << "Map Type: " << mapTypeStr(mType) << endl;
205  }
206  if ( update_basis(sock, vwl, eType, bType, bOrder, mType,
207  dType, defData, visualization) )
208  {
209  cerr << "Invalid combination of basis info (try again)" << endl;
210  }
211 
212  if (!visualization) { break; }
213 
214  print_char = false;
215  cout << endl;
216  cout << "What would you like to do?\n"
217  "q) Quit\n"
218  "c) Close Windows and Quit\n"
219  "e) Change Element Type\n"
220  "b) Change Basis Type\n";
221  if ( bType == 'h' || bType == 'p' || bType == 'n' || bType == 'r' ||
222  bType == 'l' || bType == 'f' || bType == 'g' )
223  {
224  cout << "o) Change Basis Order\n";
225  }
226  // The following is disabled pending updates to GLVis
227  if ( bType == 'l' && false )
228  {
229  cout << "m) Change Map Type\n";
230  }
231  cout << "t) Transform Element\n";
232  cout << "--> " << flush;
233  char mk;
234  cin >> mk;
235 
236  if (mk == 'q')
237  {
238  break;
239  }
240  if (mk == 'c')
241  {
242  for (unsigned int i=0; i<sock.size(); i++)
243  {
244  *sock[i] << "keys q";
245  }
246  break;
247  }
248  if (mk == 'e')
249  {
250  eInt = 0;
251  cout << "valid element types:\n";
252  if ( basisIs1D(bType) )
253  {
254  cout <<
255  "1) Segment\n";
256  }
257  if ( basisIs2D(bType) )
258  {
259  cout <<
260  "2) Triangle\n"
261  "3) Quadrilateral\n";
262  }
263  if ( basisIs3D(bType) )
264  {
265  cout <<
266  "4) Tetrahedron\n"
267  "5) Hexahedron\n";
268  }
269  cout << "enter new element type --> " << flush;
270  cin >> eInt;
271  if ( eInt <= 0 || eInt > 5 )
272  {
273  cout << "invalid element type \"" << eInt << "\"" << endl << flush;
274  }
275  else if ( (elemIs1D((Element::Type)eInt) && basisIs1D(bType)) ||
276  (elemIs2D((Element::Type)eInt) && basisIs2D(bType)) ||
277  (elemIs3D((Element::Type)eInt) && basisIs3D(bType)) )
278  {
279  if ( (elemIs1D((Element::Type)eInt) && !elemIs1D(eType)) ||
280  (elemIs2D((Element::Type)eInt) && !elemIs2D(eType)) ||
281  (elemIs3D((Element::Type)eInt) && !elemIs3D(eType)) )
282  {
283  dType = Deformation::INVALID;
284  }
285  eType = (Element::Type)eInt;
286 
287  print_char = true;
288  }
289  else
290  {
291  cout << "invalid element type \"" << eInt <<
292  "\" for basis type \"" << basisTypeStr(bType) << "\"." << endl;
293  }
294  }
295  if (mk == 'b')
296  {
297  char bChar = 0;
298  cout << "valid basis types:\n";
299  cout << "h) H1 Finite Element\n";
300  cout << "p) H1 Positive Finite Element\n";
301  if ( elemIs2D(eType) || elemIs3D(eType) )
302  {
303  cout << "n) Nedelec Finite Element\n";
304  cout << "r) Raviart-Thomas Finite Element\n";
305  }
306  cout << "l) L2 Finite Element\n";
307  if ( elemIs1D(eType) || elemIs2D(eType) )
308  {
309  cout << "c) Crouzeix-Raviart Finite Element\n";
310  }
311  cout << "f) Fixed Order Continuous Finite Element\n";
312  if ( elemIs2D(eType) )
313  {
314  cout << "g) Gauss Discontinuous Finite Element\n";
315  }
316  cout << "enter new basis type --> " << flush;
317  cin >> bChar;
318  if ( bChar == 'h' || bChar == 'p' || bChar == 'l' || bChar == 'f' ||
319  ((bChar == 'n' || bChar == 'r') &&
320  (elemIs2D(eType) || elemIs3D(eType))) ||
321  (bChar == 'c' && (elemIs1D(eType) || elemIs2D(eType))) ||
322  (bChar == 'g' && elemIs2D(eType)))
323  {
324  bType = bChar;
325  if ( bType == 'h' )
326  {
327  mType = FiniteElement::VALUE;
328  }
329  else if ( bType == 'p' )
330  {
331  mType = FiniteElement::VALUE;
332  }
333  else if ( bType == 'n' )
334  {
335  mType = FiniteElement::H_CURL;
336  }
337  else if ( bType == 'r' )
338  {
339  mType = FiniteElement::H_DIV;
340  }
341  else if ( bType == 'l' )
342  {
343  if ( mType != FiniteElement::VALUE &&
344  mType != FiniteElement::INTEGRAL )
345  {
346  mType = FiniteElement::VALUE;
347  }
348  }
349  else if ( bType == 'c' )
350  {
351  bOrder = 1;
352  mType = FiniteElement::VALUE;
353  }
354  else if ( bType == 'f' )
355  {
356  if ( bOrder < 1 || bOrder > 3)
357  {
358  bOrder = 1;
359  }
360  mType = FiniteElement::VALUE;
361  }
362  else if ( bType == 'g' )
363  {
364  if ( bOrder < 1 || bOrder > 2)
365  {
366  bOrder = 1;
367  }
368  mType = FiniteElement::VALUE;
369  }
370  print_char = true;
371  }
372  else
373  {
374  cout << "invalid basis type \"" << bChar << "\"." << endl;
375  }
376  }
377  if (mk == 'm' && bType == 'l')
378  {
379  int mInt = 0;
380  cout << "valid map types:\n"
381  "0) VALUE\n"
382  "1) INTEGRAL\n";
383  cout << "enter new map type --> " << flush;
384  cin >> mInt;
385  if (mInt >=0 && mInt <= 1)
386  {
387  mType = mInt;
388  print_char = true;
389  }
390  else
391  {
392  cout << "invalid map type \"" << mInt << "\"." << endl;
393  }
394  }
395  if (mk == 'o')
396  {
397  int oInt = 1;
398  int oMin = ( bType == 'h' || bType == 'p' || bType == 'n' ||
399  bType == 'f' || bType == 'g')?1:0;
400  int oMax = -1;
401  switch (bType)
402  {
403  case 'g':
404  oMax = 2;
405  break;
406  case 'f':
407  oMax = 3;
408  break;
409  default:
410  oMax = -1;
411  }
412  cout << "basis function order must be >= " << oMin;
413  if ( oMax >= 0 )
414  {
415  cout << " and <= " << oMax;
416  }
417  cout << endl;
418  cout << "enter new basis function order --> " << flush;
419  cin >> oInt;
420  if ( oInt >= oMin && oInt <= (oMax>=0)?oMax:oInt )
421  {
422  bOrder = oInt;
423  print_char = true;
424  }
425  else
426  {
427  cout << "invalid basis order \"" << oInt << "\"." << endl;
428  }
429  }
430  if (mk == 't')
431  {
432  cout << "transformation options:\n";
433  cout << "r) reset to reference element\n";
434  cout << "u) uniform scaling\n";
435  if ( elemIs2D(eType) || elemIs3D(eType) )
436  {
437  cout << "c) compression\n";
438  cout << "s) shear\n";
439  }
440  cout << "enter transformation type --> " << flush;
441  char tk;
442  cin >> tk;
443  if (tk == 'r')
444  {
445  dType = Deformation::INVALID;
446  }
447  else if (tk == 'u')
448  {
449  cout << "enter scaling constant --> " << flush;
450  cin >> defData.uniformScale;
451  if ( defData.uniformScale > 0.0 )
452  {
453  dType = Deformation::UNIFORM;
454  }
455  }
456  else if (tk == 'c' && !elemIs1D(eType))
457  {
458  int dim = elemIs2D(eType)?2:3;
459  cout << "enter compression factor --> " << flush;
460  cin >> defData.squeezeFactor;
461  cout << "enter compression axis (0-" << dim-1 << ") --> " << flush;
462  cin >> defData.squeezeAxis;
463 
464  if ( defData.squeezeFactor > 0.0 &&
465  (defData.squeezeAxis >= 0 && defData.squeezeAxis < dim))
466  {
467  dType = Deformation::SQUEEZE;
468  }
469  }
470  else if (tk == 's' && !elemIs1D(eType))
471  {
472  int dim = elemIs2D(eType)?2:3;
473  cout << "enter shear vector (components separated by spaces) --> "
474  << flush;
475  defData.shearVec.SetSize(dim);
476  for (int i=0; i<dim; i++)
477  {
478  cin >> defData.shearVec[i];
479  }
480  cout << "enter shear axis (0-" << dim-1 << ") --> " << flush;
481  cin >> defData.shearAxis;
482 
483  if ( defData.shearAxis >= 0 && defData.shearAxis < dim )
484  {
485  dType = Deformation::SHEAR;
486  }
487  }
488  }
489  }
490 
491  // Cleanup
492  for (unsigned int i=0; i<sock.size(); i++)
493  {
494  delete sock[i];
495  }
496 
497  // Exit
498  return 0;
499 }
500 
501 string elemTypeStr(const Element::Type & eType)
502 {
503  switch (eType)
504  {
505  case Element::POINT:
506  return "POINT";
507  case Element::SEGMENT:
508  return "SEGMENT";
509  case Element::TRIANGLE:
510  return "TRIANGLE";
511  case Element::QUADRILATERAL:
512  return "QUADRILATERAL";
513  case Element::TETRAHEDRON:
514  return "TETRAHEDRON";
515  case Element::HEXAHEDRON:
516  return "HEXAHEDRON";
517  default:
518  return "INVALID";
519  };
520 }
521 
522 bool
523 elemIs1D(const Element::Type & eType)
524 {
525  return eType == Element::SEGMENT;
526 }
527 
528 bool
529 elemIs2D(const Element::Type & eType)
530 {
531  return eType == Element::TRIANGLE || eType == Element::QUADRILATERAL;
532 }
533 
534 bool
535 elemIs3D(const Element::Type & eType)
536 {
537  return eType == Element::TETRAHEDRON || eType == Element::HEXAHEDRON;
538 }
539 
540 string
541 basisTypeStr(char bType)
542 {
543  switch (bType)
544  {
545  case 'h':
546  return "Continuous (H1)";
547  case 'p':
548  return "Continuous Positive (H1)";
549  case 'n':
550  return "Nedelec";
551  case 'r':
552  return "Raviart-Thomas";
553  case 'l':
554  return "Discontinuous (L2)";
555  case 'f':
556  return "Fixed Order Continuous";
557  case 'g':
558  return "Gaussian Discontinuous";
559  case 'c':
560  return "Crouzeix-Raviart";
561  default:
562  return "INVALID";
563  };
564 }
565 
566 bool
567 basisIs1D(char bType)
568 {
569  return bType == 'h' || bType == 'p' || bType == 'l' || bType == 'c' ||
570  bType == 'f';
571 }
572 
573 bool
574 basisIs2D(char bType)
575 {
576  return bType == 'h' || bType == 'p' || bType == 'n' || bType == 'r' ||
577  bType == 'l' || bType == 'c' || bType == 'f' || bType == 'g';
578 }
579 
580 bool
581 basisIs3D(char bType)
582 {
583  return bType == 'h' || bType == 'p' || bType == 'n' || bType == 'r' ||
584  bType == 'f' || bType == 'l';
585 }
586 
587 string
588 mapTypeStr(int mType)
589 {
590  switch (mType)
591  {
592  case FiniteElement::VALUE:
593  return "VALUE";
594  case FiniteElement::H_CURL:
595  return "H_CURL";
596  case FiniteElement::H_DIV:
597  return "H_DIV";
598  case FiniteElement::INTEGRAL:
599  return "INTEGRAL";
600  default:
601  return "INVALID";
602  }
603 }
604 
605 void
607  const IntegrationPoint &ip)
608 {
609  Vector u(dim_);
610  T.Transform(ip, u);
611 
612  switch (dim_)
613  {
614  case 1:
615  Def1D(u, v);
616  break;
617  case 2:
618  Def2D(u, v);
619  break;
620  case 3:
621  Def3D(u, v);
622  break;
623  }
624 }
625 
626 void
627 Deformation::Def1D(const Vector & u, Vector & v)
628 {
629  v = u;
630  if ( dType_ == UNIFORM )
631  {
632  v *= data_.uniformScale;
633  }
634 }
635 
636 void
637 Deformation::Def2D(const Vector & u, Vector & v)
638 {
639  switch (dType_)
640  {
641  case UNIFORM:
642  v = u;
643  v *= data_.uniformScale;
644  break;
645  case SQUEEZE:
646  v = u;
647  v[ data_.squeezeAxis ] /= data_.squeezeFactor;
648  v[(data_.squeezeAxis+1)%2] *= data_.squeezeFactor;
649  break;
650  case SHEAR:
651  v = u;
652  v.Add(v[data_.shearAxis], data_.shearVec);
653  break;
654  default:
655  v = u;
656  }
657 }
658 
659 void
660 Deformation::Def3D(const Vector & u, Vector & v)
661 {
662  switch (dType_)
663  {
664  case UNIFORM:
665  v = u;
666  v *= data_.uniformScale;
667  break;
668  case SQUEEZE:
669  v = u;
670  v[ data_.squeezeAxis ] /= data_.squeezeFactor;
671  v[(data_.squeezeAxis+1)%2] *= sqrt(data_.squeezeFactor);
672  v[(data_.squeezeAxis+2)%2] *= sqrt(data_.squeezeFactor);
673  break;
674  case SHEAR:
675  v = u;
676  v.Add(v[data_.shearAxis], data_.shearVec);
677  break;
678  default:
679  v = u;
680  }
681 }
682 
683 int
684 update_basis(vector<socketstream*> & sock, const VisWinLayout & vwl,
685  Element::Type e, char bType, int bOrder, int mType,
686  Deformation::DefType dType, const DeformationData & defData,
687  bool visualization)
688 {
689  bool vec = false;
690 
691  Mesh *mesh;
692  ElementMeshStream imesh(e);
693  if (!imesh)
694  {
695  {
696  cerr << "\nProblem with meshstream object\n" << endl;
697  }
698  return 2;
699  }
700  mesh = new Mesh(imesh, 1, 1);
701  int dim = mesh->Dimension();
702 
703  if ( dType != Deformation::INVALID )
704  {
705  Deformation defCoef(dim, dType, defData);
706  mesh->Transform(defCoef);
707  }
708 
709  FiniteElementCollection * FEC = NULL;
710  switch (bType)
711  {
712  case 'h':
713  FEC = new H1_FECollection(bOrder, dim);
714  vec = false;
715  break;
716  case 'p':
717  FEC = new H1Pos_FECollection(bOrder, dim);
718  vec = false;
719  break;
720  case 'n':
721  FEC = new ND_FECollection(bOrder, dim);
722  vec = true;
723  break;
724  case 'r':
725  FEC = new RT_FECollection(bOrder-1, dim);
726  vec = true;
727  break;
728  case 'l':
729  FEC = new L2_FECollection(bOrder, dim, BasisType::GaussLegendre,
730  mType);
731  vec = false;
732  break;
733  case 'c':
734  FEC = new CrouzeixRaviartFECollection();
735  break;
736  case 'f':
737  if ( bOrder == 1 )
738  {
739  FEC = new LinearFECollection();
740  }
741  else if ( bOrder == 2 )
742  {
743  FEC = new QuadraticFECollection();
744  }
745  else if ( bOrder == 3 )
746  {
747  FEC = new CubicFECollection();
748  }
749  break;
750  case 'g':
751  if ( bOrder == 1 )
752  {
754  }
755  else if ( bOrder == 2 )
756  {
758  }
759  break;
760  }
761  if ( FEC == NULL)
762  {
763  delete mesh;
764  return 1;
765  }
766 
767  FiniteElementSpace FESpace(mesh, FEC);
768 
769  int ndof = FESpace.GetVSize();
770 
771  Array<int> vdofs;
772  FESpace.GetElementVDofs(0,vdofs);
773 
774  char vishost[] = "localhost";
775  int visport = 19916;
776 
777  int offx = vwl.w+10, offy = vwl.h+45; // window offsets
778 
779  for (unsigned int i=0; i<sock.size(); i++)
780  {
781  *sock[i] << "keys q";
782  delete sock[i];
783  }
784 
785  sock.resize(ndof);
786  for (int i=0; i<ndof; i++)
787  {
788  sock[i] = new socketstream; sock[i]->precision(8);
789  }
790 
791  GridFunction ** x = new GridFunction*[ndof];
792  for (int i=0; i<ndof; i++)
793  {
794  x[i] = new GridFunction(&FESpace);
795  *x[i] = 0.0;
796  if ( vdofs[i] < 0 )
797  {
798  (*x[i])(-1-vdofs[i]) = -1.0;
799  }
800  else
801  {
802  (*x[i])(vdofs[i]) = 1.0;
803  }
804  }
805 
806  int ref = 0;
807  int exOrder = 0;
808  if ( bType == 'n' ) { exOrder++; }
809  if ( bType == 'r' ) { exOrder += 2; }
810  while ( 1<<ref < bOrder + exOrder || ref == 0 )
811  {
812  mesh->UniformRefinement();
813  FESpace.Update();
814 
815  for (int i=0; i<ndof; i++)
816  {
817  x[i]->Update();
818  }
819  ref++;
820  }
821 
822  for (int i=0; i<ndof; i++)
823  {
824  ostringstream oss;
825  oss << "DoF " << i + 1;
826  if (visualization)
827  {
828  VisualizeField(*sock[i], vishost, visport, *x[i], oss.str().c_str(),
829  (i % vwl.nx) * offx, ((i / vwl.nx) % vwl.ny) * offy,
830  vwl.w, vwl.h,
831  "aaAc", vec);
832  }
833  }
834 
835  for (int i=0; i<ndof; i++)
836  {
837  delete x[i];
838  }
839  delete [] x;
840 
841  delete FEC;
842  delete mesh;
843 
844  return 0;
845 }
bool elemIs3D(const Element::Type &eType)
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:661
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:347
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:1942
bool elemIs1D(const Element::Type &eType)
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:172
int offx
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:332
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:9180
int main(int argc, char *argv[])
Definition: ex1.cpp:58
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:404
string mapTypeStr(int mType)
int update_basis(vector< socketstream * > &sock, const VisWinLayout &vwl, Element::Type e, char bType, int bOrder, int mType, Deformation::DefType dType, const DeformationData &defData, bool visualization)
int dim
Definition: ex3.cpp:48
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
bool basisIs1D(char bType)
virtual void Eval(DenseMatrix &M, ElementTransformation &T, const IntegrationRule &ir)
Evaluate the vector coefficient in the element described by T at all points of ir, storing the result in M.
Definition: coefficient.cpp:90
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
string elemTypeStr(const Element::Type &eType)
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:431
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:180
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Definition: optparser.hpp:76
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:152
bool basisIs3D(char bType)
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:585
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:190
Class for integration point with weight.
Definition: intrules.hpp:25
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
bool elemIs2D(const Element::Type &eType)
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:357
bool basisIs2D(char bType)
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:240
Vector data type.
Definition: vector.hpp:48
virtual void Transform(const IntegrationPoint &, Vector &)=0
int offy
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Definition: fem_extras.cpp:94
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:134
string basisTypeStr(char bType)
bool Good() const
Definition: optparser.hpp:122