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