31 #include "../common/mfem-common.hpp"
37 using namespace mfem::common;
49 struct DeformationData
77 enum DefType {INVALID, UNIFORM, SQUEEZE, SHEAR};
79 Deformation(
int dim, DefType dType,
const DeformationData & data)
83 using VectorCoefficient::Eval;
91 const DeformationData & data_;
106 int update_basis(vector<socketstream*> & sock,
const VisWinLayout & vwl,
108 Deformation::DefType dType,
const DeformationData & defData,
109 bool visualization,
int &onlySome);
111 int main(
int argc,
char *argv[])
128 Deformation::DefType dType = Deformation::INVALID;
129 DeformationData defData;
131 bool visualization =
true;
134 vector<socketstream*> sock;
137 args.
AddOption(&eInt,
"-e",
"--elem-type",
138 "Element Type: (1-Segment, 2-Triangle, 3-Quadrilateral, "
139 "4-Tetrahedron, 5-Hexahedron)");
140 args.
AddOption(&bInt,
"-b",
"--basis-type",
141 "Basis Function Type (0-H1, 1-Nedelec, 2-Raviart-Thomas, "
142 "3-L2, 4-Fixed Order Cont.,\n\t5-Gaussian Discontinuous (2D),"
143 " 6-Crouzeix-Raviart, 7-Serendipity)");
144 args.
AddOption(&bOrder,
"-o",
"--order",
"Basis function order");
145 args.
AddOption(&vwl.nx,
"-nx",
"--num-win-x",
146 "Number of Viz windows in X");
147 args.
AddOption(&vwl.ny,
"-ny",
"--num-win-y",
148 "Number of Viz windows in y");
150 "Width of Viz windows");
152 "Height of Viz windows");
153 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
154 "--no-visualization",
155 "Enable or disable GLVis visualization.");
156 args.
AddOption(&onlySome,
"-only",
"--onlySome",
157 "Only view 10 dofs, starting with the specified one.");
167 if ( eInt > 0 && eInt < 6 )
202 bool print_char =
true;
208 cout <<
"Element Type: " <<
elemTypeStr(eType) << endl;
210 cout <<
"Basis function order: " << bOrder << endl;
211 cout <<
"Map Type: " <<
mapTypeStr(mType) << endl;
213 if (
update_basis(sock, vwl, eType, bType, bOrder, mType,
214 dType, defData, visualization, onlySome) )
216 cerr <<
"Invalid combination of basis info (try again)" << endl;
219 if (!visualization) {
break; }
223 cout <<
"What would you like to do?\n"
225 "c) Close Windows and Quit\n"
226 "e) Change Element Type\n"
227 "b) Change Basis Type\n";
228 if ( bType ==
'h' || bType ==
'p' || bType ==
'n' || bType ==
'r' ||
229 bType ==
'l' || bType ==
'f' || bType ==
'g' || bType ==
's')
231 cout <<
"o) Change Basis Order\n";
234 if ( bType ==
'l' &&
false )
236 cout <<
"m) Change Map Type\n";
238 cout <<
"t) Transform Element\n";
239 cout <<
"--> " << flush;
249 for (
unsigned int i=0; i<sock.size(); i++)
251 *sock[i] <<
"keys q";
258 cout <<
"valid element types:\n";
268 "3) Quadrilateral\n";
276 cout <<
"enter new element type --> " << flush;
278 if ( eInt <= 0 || eInt > 5 )
280 cout <<
"invalid element type \"" << eInt <<
"\"" << endl << flush;
290 dType = Deformation::INVALID;
298 cout <<
"invalid element type \"" << eInt <<
299 "\" for basis type \"" <<
basisTypeStr(bType) <<
"\"." << endl;
305 cout <<
"valid basis types:\n";
306 cout <<
"h) H1 Finite Element\n";
307 cout <<
"p) H1 Positive Finite Element\n";
310 cout <<
"s) H1 Serendipity Finite Element\n";
311 cout <<
"n) Nedelec Finite Element\n";
312 cout <<
"r) Raviart-Thomas Finite Element\n";
314 cout <<
"l) L2 Finite Element\n";
317 cout <<
"c) Crouzeix-Raviart Finite Element\n";
319 cout <<
"f) Fixed Order Continuous Finite Element\n";
322 cout <<
"g) Gauss Discontinuous Finite Element\n";
324 cout <<
"enter new basis type --> " << flush;
326 if (bChar ==
'h' || bChar ==
'p' || bChar ==
'l' || bChar ==
'f' ||
335 mType = FiniteElement::VALUE;
337 else if ( bType ==
'p' )
339 mType = FiniteElement::VALUE;
341 else if (bType ==
's')
343 mType = FiniteElement::VALUE;
345 else if ( bType ==
'n' )
347 mType = FiniteElement::H_CURL;
349 else if ( bType ==
'r' )
351 mType = FiniteElement::H_DIV;
353 else if ( bType ==
'l' )
355 if ( mType != FiniteElement::VALUE &&
356 mType != FiniteElement::INTEGRAL )
358 mType = FiniteElement::VALUE;
361 else if ( bType ==
'c' )
364 mType = FiniteElement::VALUE;
366 else if ( bType ==
'f' )
368 if ( bOrder < 1 || bOrder > 3)
372 mType = FiniteElement::VALUE;
374 else if ( bType ==
'g' )
376 if ( bOrder < 1 || bOrder > 2)
380 mType = FiniteElement::VALUE;
386 cout <<
"invalid basis type \"" << bChar <<
"\"." << endl;
389 if (mk ==
'm' && bType ==
'l')
392 cout <<
"valid map types:\n"
395 cout <<
"enter new map type --> " << flush;
397 if (mInt >=0 && mInt <= 1)
404 cout <<
"invalid map type \"" << mInt <<
"\"." << endl;
410 int oMin = ( bType ==
'h' || bType ==
'p' || bType ==
'n' ||
411 bType ==
'f' || bType ==
'g' || bType ==
's')?1:0;
424 cout <<
"basis function order must be >= " << oMin;
427 cout <<
" and <= " << oMax;
430 cout <<
"enter new basis function order --> " << flush;
432 if ( oInt >= oMin && oInt <= (oMax>=0)?oMax:oInt )
439 cout <<
"invalid basis order \"" << oInt <<
"\"." << endl;
444 cout <<
"transformation options:\n";
445 cout <<
"r) reset to reference element\n";
446 cout <<
"u) uniform scaling\n";
449 cout <<
"c) compression\n";
450 cout <<
"s) shear\n";
452 cout <<
"enter transformation type --> " << flush;
457 dType = Deformation::INVALID;
461 cout <<
"enter scaling constant --> " << flush;
462 cin >> defData.uniformScale;
463 if ( defData.uniformScale > 0.0 )
465 dType = Deformation::UNIFORM;
468 else if (tk ==
'c' && !
elemIs1D(eType))
471 cout <<
"enter compression factor --> " << flush;
472 cin >> defData.squeezeFactor;
473 cout <<
"enter compression axis (0-" << dim-1 <<
") --> " << flush;
474 cin >> defData.squeezeAxis;
476 if ( defData.squeezeFactor > 0.0 &&
477 (defData.squeezeAxis >= 0 && defData.squeezeAxis < dim))
479 dType = Deformation::SQUEEZE;
482 else if (tk ==
's' && !
elemIs1D(eType))
485 cout <<
"enter shear vector (components separated by spaces) --> "
487 defData.shearVec.SetSize(dim);
488 for (
int i=0; i<
dim; i++)
490 cin >> defData.shearVec[i];
492 cout <<
"enter shear axis (0-" << dim-1 <<
") --> " << flush;
493 cin >> defData.shearAxis;
495 if ( defData.shearAxis >= 0 && defData.shearAxis < dim )
497 dType = Deformation::SHEAR;
504 for (
unsigned int i=0; i<sock.size(); i++)
519 case Element::SEGMENT:
521 case Element::TRIANGLE:
523 case Element::QUADRILATERAL:
524 return "QUADRILATERAL";
525 case Element::TETRAHEDRON:
526 return "TETRAHEDRON";
527 case Element::HEXAHEDRON:
537 return eType == Element::SEGMENT;
543 return eType == Element::TRIANGLE || eType == Element::QUADRILATERAL;
549 return eType == Element::TETRAHEDRON || eType == Element::HEXAHEDRON;
558 return "Continuous (H1)";
560 return "Continuous Positive (H1)";
562 return "Continuous Serendipity (H1)";
566 return "Raviart-Thomas";
568 return "Discontinuous (L2)";
570 return "Fixed Order Continuous";
572 return "Gaussian Discontinuous";
574 return "Crouzeix-Raviart";
583 return bType ==
'h' || bType ==
'p' || bType ==
'l' || bType ==
'c' ||
590 return bType ==
'h' || bType ==
'p' || bType ==
'n' || bType ==
'r' ||
591 bType ==
'l' || bType ==
'c' || bType ==
'f' || bType ==
'g' ||
598 return bType ==
'h' || bType ==
'p' || bType ==
'n' || bType ==
'r' ||
599 bType ==
'f' || bType ==
'l';
607 case FiniteElement::VALUE:
609 case FiniteElement::H_CURL:
611 case FiniteElement::H_DIV:
613 case FiniteElement::INTEGRAL:
645 if ( dType_ == UNIFORM )
647 v *= data_.uniformScale;
658 v *= data_.uniformScale;
662 v[ data_.squeezeAxis ] /= data_.squeezeFactor;
663 v[(data_.squeezeAxis+1)%2] *= data_.squeezeFactor;
667 v.
Add(v[data_.shearAxis], data_.shearVec);
681 v *= data_.uniformScale;
685 v[ data_.squeezeAxis ] /= data_.squeezeFactor;
686 v[(data_.squeezeAxis+1)%2] *= sqrt(data_.squeezeFactor);
687 v[(data_.squeezeAxis+2)%2] *= sqrt(data_.squeezeFactor);
691 v.
Add(v[data_.shearAxis], data_.shearVec);
701 Deformation::DefType dType,
const DeformationData & defData,
702 bool visualization,
int &onlySome)
711 cerr <<
"\nProblem with meshstream object\n" << endl;
715 mesh =
new Mesh(imesh, 1, 1);
718 if ( dType != Deformation::INVALID )
720 Deformation defCoef(dim, dType, defData);
767 else if ( bOrder == 2 )
771 else if ( bOrder == 3 )
781 else if ( bOrder == 2 )
803 int offx = vwl.w+10,
offy = vwl.h+45;
805 for (
unsigned int i=0; i<sock.size(); i++)
807 *sock[i] <<
"keys q";
812 for (
int i=0; i<ndof; i++)
818 for (
int i=0; i<ndof; i++)
824 (*x[i])(-1-vdofs[i]) = -1.0;
828 (*x[i])(vdofs[i]) = 1.0;
834 if ( bType ==
'n' ) { exOrder++; }
835 if ( bType ==
'r' ) { exOrder += 2; }
836 while ( 1<<ref < bOrder + exOrder || ref == 0 )
841 for (
int i=0; i<ndof; i++)
849 if (ndof > 25 && onlySome == -1)
852 cout <<
"There are more than 25 windows to open.\n"
853 <<
"Only showing Dofs 1-10 to avoid crashing.\n"
854 <<
"Use the option -only N to show Dofs N to N+9 instead.\n";
857 for (
int i = 0; i < stopAt; i++)
859 if (i ==0 && onlySome > 0 && onlySome <ndof)
862 stopAt = min(ndof,onlySome+9);
866 oss <<
"DoF " << i + 1;
869 VisualizeField(*sock[i], vishost, visport, *x[i], oss.str().c_str(),
870 (i % vwl.nx) * offx, ((i / vwl.nx) % vwl.ny) *
offy,
876 for (
int i=0; i<ndof; i++)
bool elemIs3D(const Element::Type &eType)
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Class for grid function - Vector with associated FE space.
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Base class for vector Coefficients that optionally depend on time and space.
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'th element.
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)
Piecewise-(bi/tri)linear continuous finite elements.
void Transform(void(*f)(const Vector &, Vector &))
int main(int argc, char *argv[])
Piecewise-(bi)cubic continuous finite elements.
string mapTypeStr(int mType)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
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.
Type
Constants for the classes derived from Element.
string elemTypeStr(const Element::Type &eType)
void PrintUsage(std::ostream &out) const
Print the usage message.
Crouzeix-Raviart nonconforming elements in 2D.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
bool basisIs3D(char bType)
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
Print the options.
bool elemIs2D(const Element::Type &eType)
Piecewise-(bi)quadratic continuous finite elements.
bool basisIs2D(char bType)
Arbitrary order H(curl)-conforming Nedelec finite elements.
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)
Arbitrary order H1-conforming (continuous) finite elements.
Arbitrary order "L2-conforming" discontinuous finite elements.
string basisTypeStr(char bType)
bool Good() const
Return true if the command line options were parsed successfully.