MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nurbs.hpp
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 #ifndef MFEM_NURBS
13 #define MFEM_NURBS
14 
15 #include "../config/config.hpp"
16 #include "../general/table.hpp"
17 #include "../linalg/vector.hpp"
18 #include "element.hpp"
19 #include "mesh.hpp"
20 #ifdef MFEM_USE_MPI
21 #include "../general/communication.hpp"
22 #endif
23 #include <iostream>
24 
25 namespace mfem
26 {
27 
28 class GridFunction;
29 
30 
32 {
33 protected:
34  static const int MaxOrder;
35 
38 
39 public:
40  /// Create KnotVector
41  KnotVector() { }
42  KnotVector(std::istream &input);
43  KnotVector(int Order_, int NCP);
44  KnotVector(const KnotVector &kv) { (*this) = kv; }
45 
46  KnotVector &operator=(const KnotVector &kv);
47 
48  int GetNE() const { return NumOfElements; }
49  int GetNKS() const { return NumOfControlPoints - Order; }
50  int GetNCP() const { return NumOfControlPoints; }
51  int GetOrder() const { return Order; }
52  int Size() const { return knot.Size(); }
53 
54  /// Count the number of elements
55  void GetElements();
56 
57  bool isElement(int i) const { return (knot(Order+i) != knot(Order+i+1)); }
58 
59  double getKnotLocation(double xi, int ni) const
60  { return (xi*knot(ni+1) + (1. - xi)*knot(ni)); }
61 
62  int findKnotSpan(double u) const;
63 
64  void CalcShape (Vector &shape, int i, double xi) const;
65  void CalcDShape(Vector &grad, int i, double xi) const;
66 
67  void Difference(const KnotVector &kv, Vector &diff) const;
68  void UniformRefinement(Vector &newknots) const;
69  /** Return a new KnotVector with elevated degree by repeating the endpoints
70  of the knot vector. */
71  KnotVector *DegreeElevate(int t) const;
72 
73  void Flip();
74 
75  void Print(std::ostream &out) const;
76 
77  /// Destroys KnotVector
79 
80  double &operator[](int i) { return knot(i); }
81  const double &operator[](int i) const { return knot(i); }
82 };
83 
84 
86 {
87 protected:
88  int ni, nj, nk, Dim;
89  double *data; // the layout of data is: (Dim x ni x nj x nk)
90 
92 
93  int sd, nd;
94 
95  void swap(NURBSPatch *np);
96 
97  // Special B-NET access functions
98  int SetLoopDirection(int dir);
99  inline double &operator()(int i, int j);
100  inline const double &operator()(int i, int j) const;
101 
102  void init(int dim_);
103 
104  NURBSPatch(NURBSPatch *parent, int dir, int Order, int NCP);
105 
106 public:
107  NURBSPatch(const NURBSPatch &orig);
108  NURBSPatch(std::istream &input);
109  NURBSPatch(const KnotVector *kv0, const KnotVector *kv1, int dim_);
110  NURBSPatch(const KnotVector *kv0, const KnotVector *kv1,
111  const KnotVector *kv2, int dim_);
113 
114  ~NURBSPatch();
115 
116  void Print(std::ostream &out) const;
117 
118  void DegreeElevate(int dir, int t);
119  void KnotInsert (int dir, const KnotVector &knot);
120  void KnotInsert (int dir, const Vector &knot);
121 
122  void KnotInsert(Array<KnotVector *> &knot);
123  void DegreeElevate(int t);
124  void UniformRefinement();
125 
126  // Return the number of components stored in the NURBSPatch
127  int GetNC() const { return Dim; }
128  int GetNKV() const { return kv.Size(); }
129  KnotVector *GetKV(int i) { return kv[i]; }
130 
131  // Standard B-NET access functions
132  inline double &operator()(int i, int j, int l);
133  inline const double &operator()(int i, int j, int l) const;
134 
135  inline double &operator()(int i, int j, int k, int l);
136  inline const double &operator()(int i, int j, int k, int l) const;
137 
138  static void Get3DRotationMatrix(double n[], double angle, double r,
139  DenseMatrix &T);
140  void FlipDirection(int dir);
141  void SwapDirections(int dir1, int dir2);
142  void Rotate3D(double normal[], double angle);
143  int MakeUniformDegree(int degree = -1);
144  friend NURBSPatch *Interpolate(NURBSPatch &p1, NURBSPatch &p2);
145  friend NURBSPatch *Revolve3D(NURBSPatch &patch, double n[], double ang,
146  int times);
147 };
148 
149 
150 #ifdef MFEM_USE_MPI
151 class ParNURBSExtension;
152 #endif
153 
154 class NURBSPatchMap;
155 
156 
158 {
159 #ifdef MFEM_USE_MPI
160  friend class ParNURBSExtension;
161 #endif
162  friend class NURBSPatchMap;
163 
164 protected:
165  int mOrder; // see GetOrder() for description
168  // global entity counts
170  // local entity counts
173 
174  Array<int> activeVert; // activeVert[glob_vert] = loc_vert or -1
177  Array<int> activeDof; // activeDof[glob_dof] = loc_dof + 1 or 0
178 
180  int own_topo;
184 
185  // global offsets, meshOffsets == meshVertexOffsets
190 
191  // global offsets, spaceOffsets == dofOffsets
196 
198 
201  Array2D<int> el_to_IJK; // IJK are "knot-span" indices!
202  Array2D<int> bel_to_IJK; // they are NOT element indices!
203 
205 
206  inline int KnotInd(int edge) const;
207  inline KnotVector *KnotVec(int edge);
208  inline const KnotVector *KnotVec(int edge) const;
209  inline const KnotVector *KnotVec(int edge, int oedge, int *okv) const;
210 
211  void CheckPatches();
212  void CheckBdrPatches();
213 
214  void GetPatchKnotVectors (int p, Array<KnotVector *> &kv);
215  void GetPatchKnotVectors (int p, Array<const KnotVector *> &kv) const;
217  void GetBdrPatchKnotVectors(int p, Array<const KnotVector *> &kv) const;
218 
219  void SetOrderFromOrders();
221 
222  // also count the global NumOfVertices and the global NumOfDofs
223  void GenerateOffsets();
224  // count the global NumOfElements
225  void CountElements();
226  // count the global NumOfBdrElements
227  void CountBdrElements();
228 
229  // generate the mesh elements
230  void Get2DElementTopo(Array<Element *> &elements) const;
231  void Get3DElementTopo(Array<Element *> &elements) const;
232 
233  // generate the boundary mesh elements
234  void Get2DBdrElementTopo(Array<Element *> &boundary) const;
235  void Get3DBdrElementTopo(Array<Element *> &boundary) const;
236 
237 
238  // FE space generation functions
239 
240  // based on activeElem, count NumOfActiveDofs, generate el_dof,
241  // el_to_patch, el_to_IJK, activeDof map (global-to-local)
243 
244  // generate elem_to_global-dof table for the active elements
245  // define el_to_patch, el_to_IJK, activeDof (as bool)
248 
249  // call after GenerateElementDofTable
251 
252  // generate the bdr-elem_to_global-dof table for the active bdr. elements
253  // define bel_to_patch, bel_to_IJK
256 
257  // FE --> Patch translation functions
258  void GetPatchNets (const Vector &Nodes, int vdim);
259  void Get2DPatchNets(const Vector &Nodes, int vdim);
260  void Get3DPatchNets(const Vector &Nodes, int vdim);
261 
262  // Patch --> FE translation functions
263  // Side effects: delete the patches, update the weights from the patches
264  void SetSolutionVector (Vector &Nodes, int vdim);
265  void Set2DSolutionVector(Vector &Nodes, int vdim);
266  void Set3DSolutionVector(Vector &Nodes, int vdim);
267 
268  // determine activeVert, NumOfActiveVertices from the activeElem array
269  void GenerateActiveVertices();
270 
271  // determine activeBdrElem, NumOfActiveBdrElems
272  void GenerateActiveBdrElems();
273 
274  void MergeWeights(Mesh *mesh_array[], int num_pieces);
275 
276  // to be used by ParNURBSExtension constructor(s)
278 
279 public:
280  /// Copy constructor: deep copy
281  NURBSExtension(const NURBSExtension &orig);
282  /// Read-in a NURBSExtension
283  NURBSExtension(std::istream &input);
284  /** @brief Create a NURBSExtension with elevated order by repeating the
285  endpoints of the knot vectors and using uniform weights of 1. */
286  /** If a knot vector in @a parent already has order greater than or equal to
287  @a newOrder, it will be used unmodified. */
288  NURBSExtension(NURBSExtension *parent, int newOrder);
289  /** @brief Create a NURBSExtension with elevated knot vector orders (by
290  repeating the endpoints of the knot vectors and using uniform weights of
291  1) as given by the array @a newOrders. */
292  /** If a knot vector in @a parent already has order greater than or equal to
293  the corresponding entry in @a newOrder, it will be used unmodified. */
294  NURBSExtension(NURBSExtension *parent, const Array<int> &newOrders);
295  /// Construct a NURBSExtension by merging a partitioned NURBS mesh
296  NURBSExtension(Mesh *mesh_array[], int num_pieces);
297 
298  void MergeGridFunctions(GridFunction *gf_array[], int num_pieces,
299  GridFunction &merged);
300 
301  /// Destroy a NURBSExtension
302  virtual ~NURBSExtension();
303 
304  // Print functions
305  void Print(std::ostream &out) const;
306  void PrintCharacteristics(std::ostream &out) const;
307 
308  // Meta data functions
309  int Dimension() const { return patchTopo->Dimension(); }
310  int GetNP() const { return patchTopo->GetNE(); }
311  int GetNBP() const { return patchTopo->GetNBE(); }
312 
313  /// Read-only access to the orders of all knot vectors.
314  const Array<int> &GetOrders() const { return mOrders; }
315  /** @brief If all orders are identical, return that number. Otherwise, return
316  NURBSFECollection::VariableOrder. */
317  int GetOrder() const { return mOrder; }
318 
319  int GetNKV() const { return NumOfKnotVectors; }
320 
321  int GetGNV() const { return NumOfVertices; }
322  int GetNV() const { return NumOfActiveVertices; }
323  int GetGNE() const { return NumOfElements; }
324  int GetNE() const { return NumOfActiveElems; }
325  int GetGNBE() const { return NumOfBdrElements; }
326  int GetNBE() const { return NumOfActiveBdrElems; }
327 
328  int GetNTotalDof() const { return NumOfDofs; }
329  int GetNDof() const { return NumOfActiveDofs; }
330 
331  // Knotvector read-only access function
332  const KnotVector *GetKnotVector(int i) const { return knotVectors[i]; }
333 
334  // Mesh generation functions
335  void GetElementTopo (Array<Element *> &elements) const;
336  void GetBdrElementTopo(Array<Element *> &boundary) const;
337 
338  bool HavePatches() const { return (patches.Size() != 0); }
339 
342 
343  void GetVertexLocalToGlobal(Array<int> &lvert_vert);
344  void GetElementLocalToGlobal(Array<int> &lelem_elem);
345 
346  // Load functions
347  void LoadFE(int i, const FiniteElement *FE) const;
348  void LoadBE(int i, const FiniteElement *BE) const;
349 
350  const Vector &GetWeights() const { return weights; }
351  Vector &GetWeights() { return weights; }
352 
353  // Translation functions: from FE coordinates to IJK patch
354  // format and vice versa
355  void ConvertToPatches(const Vector &Nodes);
356  void SetKnotsFromPatches();
357  void SetCoordsFromPatches(Vector &Nodes);
358 
359  // Read a GridFunction written patch-by-patch, e.g. with PrintSolution().
360  void LoadSolution(std::istream &input, GridFunction &sol) const;
361  // Write a GridFunction patch-by-patch.
362  void PrintSolution(const GridFunction &sol, std::ostream &out) const;
363 
364  // Refinement methods
365  // new_degree = max(old_degree, min(old_degree + rel_degree, degree))
366  void DegreeElevate(int rel_degree, int degree = 16);
367  void UniformRefinement();
369 };
370 
371 
372 #ifdef MFEM_USE_MPI
374 {
375 private:
376  int *partitioning;
377 
378  Table *GetGlobalElementDofTable();
379  Table *Get2DGlobalElementDofTable();
380  Table *Get3DGlobalElementDofTable();
381 
382  void SetActive(const int *partitioning, const Array<bool> &active_bel);
383  void BuildGroups(const int *partitioning, const Table &elem_dof);
384 
385 public:
387 
389 
391 
392  ParNURBSExtension(MPI_Comm comm, NURBSExtension *parent, int *partitioning,
393  const Array<bool> &active_bel);
394 
395  // Create a parallel version of 'parent' with partitioning as in
396  // 'par_parent'; the 'parent' object is destroyed.
397  // The 'parent' can be either a local NURBSExtension or a global one.
399  const ParNURBSExtension *par_parent);
400 
401  virtual ~ParNURBSExtension() { delete [] partitioning; }
402 };
403 #endif
404 
405 
407 {
408 private:
409  const NURBSExtension *Ext;
410 
411  int I, J, K, pOffset, opatch;
412  Array<int> verts, edges, faces, oedge, oface;
413 
414  inline static int F(const int n, const int N)
415  { return (n < 0) ? 0 : ((n >= N) ? 2 : 1); }
416 
417  inline static int Or1D(const int n, const int N, const int Or)
418  { return (Or > 0) ? n : (N - 1 - n); }
419 
420  inline static int Or2D(const int n1, const int n2,
421  const int N1, const int N2, const int Or);
422 
423  // also set verts, edges, faces, orientations etc
424  void GetPatchKnotVectors (int p, const KnotVector *kv[]);
425  void GetBdrPatchKnotVectors(int p, const KnotVector *kv[], int *okv);
426 
427 public:
428  NURBSPatchMap(const NURBSExtension *ext) { Ext = ext; }
429 
430  int nx() { return I + 1; }
431  int ny() { return J + 1; }
432  int nz() { return K + 1; }
433 
434  void SetPatchVertexMap(int p, const KnotVector *kv[]);
435  void SetPatchDofMap (int p, const KnotVector *kv[]);
436 
437  void SetBdrPatchVertexMap(int p, const KnotVector *kv[], int *okv);
438  void SetBdrPatchDofMap (int p, const KnotVector *kv[], int *okv);
439 
440  inline int operator()(const int i) const;
441  inline int operator[](const int i) const { return (*this)(i); }
442 
443  inline int operator()(const int i, const int j) const;
444 
445  inline int operator()(const int i, const int j, const int k) const;
446 };
447 
448 
449 // Inline function implementations
450 
451 inline double &NURBSPatch::operator()(int i, int j)
452 {
453  return data[j%sd + sd*(i + (j/sd)*nd)];
454 }
455 
456 inline const double &NURBSPatch::operator()(int i, int j) const
457 {
458  return data[j%sd + sd*(i + (j/sd)*nd)];
459 }
460 
461 inline double &NURBSPatch::operator()(int i, int j, int l)
462 {
463 #ifdef MFEM_DEBUG
464  if (data == 0 || i < 0 || i >= ni || j < 0 || j >= nj || nk > 0 ||
465  l < 0 || l >= Dim)
466  {
467  mfem_error("NURBSPatch::operator() 2D");
468  }
469 #endif
470 
471  return data[(i+j*ni)*Dim+l];
472 }
473 
474 inline const double &NURBSPatch::operator()(int i, int j, int l) const
475 {
476 #ifdef MFEM_DEBUG
477  if (data == 0 || i < 0 || i >= ni || j < 0 || j >= nj || nk > 0 ||
478  l < 0 || l >= Dim)
479  {
480  mfem_error("NURBSPatch::operator() const 2D");
481  }
482 #endif
483 
484  return data[(i+j*ni)*Dim+l];
485 }
486 
487 inline double &NURBSPatch::operator()(int i, int j, int k, int l)
488 {
489 #ifdef MFEM_DEBUG
490  if (data == 0 || i < 0 || i >= ni || j < 0 || j >= nj || k < 0 ||
491  k >= nk || l < 0 || l >= Dim)
492  {
493  mfem_error("NURBSPatch::operator() 3D");
494  }
495 #endif
496 
497  return data[(i+(j+k*nj)*ni)*Dim+l];
498 }
499 
500 inline const double &NURBSPatch::operator()(int i, int j, int k, int l) const
501 {
502 #ifdef MFEM_DEBUG
503  if (data == 0 || i < 0 || i >= ni || j < 0 || j >= nj || k < 0 ||
504  k >= nk || l < 0 || l >= Dim)
505  {
506  mfem_error("NURBSPatch::operator() const 3D");
507  }
508 #endif
509 
510  return data[(i+(j+k*nj)*ni)*Dim+l];
511 }
512 
513 
514 inline int NURBSExtension::KnotInd(int edge) const
515 {
516  int kv = edge_to_knot[edge];
517  return (kv >= 0) ? kv : (-1-kv);
518 }
519 
521 {
522  return knotVectors[KnotInd(edge)];
523 }
524 
525 inline const KnotVector *NURBSExtension::KnotVec(int edge) const
526 {
527  return knotVectors[KnotInd(edge)];
528 }
529 
530 inline const KnotVector *NURBSExtension::KnotVec(int edge, int oedge, int *okv)
531 const
532 {
533  int kv = edge_to_knot[edge];
534  if (kv >= 0)
535  {
536  *okv = oedge;
537  return knotVectors[kv];
538  }
539  else
540  {
541  *okv = -oedge;
542  return knotVectors[-1-kv];
543  }
544 }
545 
546 
547 // static method
548 inline int NURBSPatchMap::Or2D(const int n1, const int n2,
549  const int N1, const int N2, const int Or)
550 {
551  // Needs testing
552  switch (Or)
553  {
554  case 0: return n1 + n2*N1;
555  case 1: return n2 + n1*N2;
556  case 2: return n2 + (N1 - 1 - n1)*N2;
557  case 3: return (N1 - 1 - n1) + n2*N1;
558  case 4: return (N1 - 1 - n1) + (N2 - 1 - n2)*N1;
559  case 5: return (N2 - 1 - n2) + (N1 - 1 - n1)*N2;
560  case 6: return (N2 - 1 - n2) + n1*N2;
561  case 7: return n1 + (N2 - 1 - n2)*N1;
562  }
563 #ifdef MFEM_DEBUG
564  mfem_error("NURBSPatchMap::Or2D");
565 #endif
566  return -1;
567 }
568 
569 inline int NURBSPatchMap::operator()(const int i) const
570 {
571  const int i1 = i - 1;
572  switch (F(i1, I))
573  {
574  case 0: return verts[0];
575  case 1: return pOffset + Or1D(i1, I, opatch);
576  case 2: return verts[1];
577  }
578 #ifdef MFEM_DEBUG
579  mfem_error("NURBSPatchMap::operator() const 1D");
580 #endif
581  return -1;
582 }
583 
584 inline int NURBSPatchMap::operator()(const int i, const int j) const
585 {
586  const int i1 = i - 1, j1 = j - 1;
587  switch (3*F(j1, J) + F(i1, I))
588  {
589  case 0: return verts[0];
590  case 1: return edges[0] + Or1D(i1, I, oedge[0]);
591  case 2: return verts[1];
592  case 3: return edges[3] + Or1D(j1, J, -oedge[3]);
593  case 4: return pOffset + Or2D(i1, j1, I, J, opatch);
594  case 5: return edges[1] + Or1D(j1, J, oedge[1]);
595  case 6: return verts[3];
596  case 7: return edges[2] + Or1D(i1, I, -oedge[2]);
597  case 8: return verts[2];
598  }
599 #ifdef MFEM_DEBUG
600  mfem_error("NURBSPatchMap::operator() const 2D");
601 #endif
602  return -1;
603 }
604 
605 inline int NURBSPatchMap::operator()(const int i, const int j, const int k)
606 const
607 {
608  // Needs testing
609  const int i1 = i - 1, j1 = j - 1, k1 = k - 1;
610  switch (3*(3*F(k1, K) + F(j1, J)) + F(i1, I))
611  {
612  case 0: return verts[0];
613  case 1: return edges[0] + Or1D(i1, I, oedge[0]);
614  case 2: return verts[1];
615  case 3: return edges[3] + Or1D(j1, J, oedge[3]);
616  case 4: return faces[0] + Or2D(i1, J - 1 - j1, I, J, oface[0]);
617  case 5: return edges[1] + Or1D(j1, J, oedge[1]);
618  case 6: return verts[3];
619  case 7: return edges[2] + Or1D(i1, I, oedge[2]);
620  case 8: return verts[2];
621  case 9: return edges[8] + Or1D(k1, K, oedge[8]);
622  case 10: return faces[1] + Or2D(i1, k1, I, K, oface[1]);
623  case 11: return edges[9] + Or1D(k1, K, oedge[9]);
624  case 12: return faces[4] + Or2D(J - 1 - j1, k1, J, K, oface[4]);
625  case 13: return pOffset + I*(J*k1 + j1) + i1;
626  case 14: return faces[2] + Or2D(j1, k1, J, K, oface[2]);
627  case 15: return edges[11] + Or1D(k1, K, oedge[11]);
628  case 16: return faces[3] + Or2D(I - 1 - i1, k1, I, K, oface[3]);
629  case 17: return edges[10] + Or1D(k1, K, oedge[10]);
630  case 18: return verts[4];
631  case 19: return edges[4] + Or1D(i1, I, oedge[4]);
632  case 20: return verts[5];
633  case 21: return edges[7] + Or1D(j1, J, oedge[7]);
634  case 22: return faces[5] + Or2D(i1, j1, I, J, oface[5]);
635  case 23: return edges[5] + Or1D(j1, J, oedge[5]);
636  case 24: return verts[7];
637  case 25: return edges[6] + Or1D(i1, I, oedge[6]);
638  case 26: return verts[6];
639  }
640 #ifdef MFEM_DEBUG
641  mfem_error("NURBSPatchMap::operator() const 3D");
642 #endif
643  return -1;
644 }
645 
646 }
647 
648 #endif
Array< KnotVector * > knotVectors
Definition: nurbs.hpp:182
Abstract class for Finite Elements.
Definition: fe.hpp:140
Array< KnotVector * > kv
Definition: nurbs.hpp:91
Array< int > f_meshOffsets
Definition: nurbs.hpp:188
KnotVector * DegreeElevate(int t) const
Definition: nurbs.cpp:57
NURBSPatch(NURBSPatch *parent, int dir, int Order, int NCP)
Definition: nurbs.cpp:404
Array2D< int > bel_to_IJK
Definition: nurbs.hpp:202
void SetCoordsFromPatches(Vector &Nodes)
Definition: nurbs.cpp:2562
int Size() const
Definition: nurbs.hpp:52
int GetGNE() const
Definition: nurbs.hpp:323
Array< int > activeVert
Definition: nurbs.hpp:174
KnotVector * KnotVec(int edge)
Definition: nurbs.hpp:520
const Vector & GetWeights() const
Definition: nurbs.hpp:350
void swap(NURBSPatch *np)
Definition: nurbs.cpp:419
int GetNV() const
Definition: nurbs.hpp:322
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void DegreeElevate(int dir, int t)
Definition: nurbs.cpp:697
int NumOfElements
Definition: nurbs.hpp:37
void LoadBE(int i, const FiniteElement *BE) const
Definition: nurbs.cpp:2530
void Generate2DElementDofTable()
Definition: nurbs.cpp:2249
int KnotInd(int edge) const
Definition: nurbs.hpp:514
void MergeGridFunctions(GridFunction *gf_array[], int num_pieces, GridFunction &merged)
Definition: nurbs.cpp:1711
Array< int > ldof_group
Definition: nurbs.hpp:388
virtual ~ParNURBSExtension()
Definition: nurbs.hpp:401
void Generate3DElementDofTable()
Definition: nurbs.cpp:2303
void KnotInsert(int dir, const KnotVector &knot)
Definition: nurbs.cpp:559
void MergeWeights(Mesh *mesh_array[], int num_pieces)
Definition: nurbs.cpp:1686
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:621
void Get3DBdrElementTopo(Array< Element * > &boundary) const
Definition: nurbs.cpp:2181
NURBSPatchMap(const NURBSExtension *ext)
Definition: nurbs.hpp:428
void UniformRefinement()
Definition: nurbs.cpp:2696
int GetNE() const
Definition: nurbs.hpp:48
void Rotate3D(double normal[], double angle)
Definition: nurbs.cpp:1024
Array< int > e_meshOffsets
Definition: nurbs.hpp:187
int GetNKV() const
Definition: nurbs.hpp:319
int GetNKS() const
Definition: nurbs.hpp:49
int operator()(const int i) const
Definition: nurbs.hpp:569
void GetElements()
Count the number of elements.
Definition: nurbs.cpp:100
void SetOrdersFromKnotVectors()
Definition: nurbs.cpp:1904
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:120
void Generate3DBdrElementDofTable()
Definition: nurbs.cpp:2431
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
int NumOfControlPoints
Definition: nurbs.hpp:37
Array< int > e_spaceOffsets
Definition: nurbs.hpp:193
void GenerateBdrElementDofTable()
Definition: nurbs.cpp:2369
const KnotVector * GetKnotVector(int i) const
Definition: nurbs.hpp:332
void GetElementLocalToGlobal(Array< int > &lelem_elem)
Definition: nurbs.cpp:2499
void init(int dim_)
Definition: nurbs.cpp:277
Array< int > edge_to_knot
Definition: nurbs.hpp:181
void GenerateActiveVertices()
Definition: nurbs.cpp:1597
Array< bool > activeElem
Definition: nurbs.hpp:175
void GetVertexLocalToGlobal(Array< int > &lvert_vert)
Definition: nurbs.cpp:2489
void Set2DSolutionVector(Vector &Nodes, int vdim)
Definition: nurbs.cpp:2812
Vector knot
Definition: nurbs.hpp:36
int GetNTotalDof() const
Definition: nurbs.hpp:328
friend NURBSPatch * Revolve3D(NURBSPatch &patch, double n[], double ang, int times)
Definition: nurbs.cpp:1120
int GetNDof() const
Definition: nurbs.hpp:329
int Dimension() const
Definition: nurbs.hpp:309
Array< int > el_to_patch
Definition: nurbs.hpp:199
void SwapDirections(int dir1, int dir2)
Definition: nurbs.cpp:960
KnotVector()
Create KnotVector.
Definition: nurbs.hpp:41
static const int MaxOrder
Definition: nurbs.hpp:34
Array< int > v_spaceOffsets
Definition: nurbs.hpp:192
void UniformRefinement(Vector &newknots) const
Definition: nurbs.cpp:86
void Generate2DBdrElementDofTable()
Definition: nurbs.cpp:2388
ParNURBSExtension(const ParNURBSExtension &orig)
Definition: nurbs.cpp:2873
void Set3DSolutionVector(Vector &Nodes, int vdim)
Definition: nurbs.cpp:2840
void GetPatchKnotVectors(int p, Array< KnotVector * > &kv)
Definition: nurbs.cpp:1810
void GenerateActiveBdrElems()
Definition: nurbs.cpp:1664
double * data
Definition: nurbs.hpp:89
void PrintSolution(const GridFunction &sol, std::ostream &out) const
Definition: nurbs.cpp:2643
void Get3DPatchNets(const Vector &Nodes, int vdim)
Definition: nurbs.cpp:2770
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:146
void Get2DPatchNets(const Vector &Nodes, int vdim)
Definition: nurbs.cpp:2743
Vector & GetWeights()
Definition: nurbs.hpp:351
bool HavePatches() const
Definition: nurbs.hpp:338
virtual ~NURBSExtension()
Destroy a NURBSExtension.
Definition: nurbs.cpp:1508
void GenerateOffsets()
Definition: nurbs.cpp:1914
Array< int > activeDof
Definition: nurbs.hpp:177
friend NURBSPatch * Interpolate(NURBSPatch &p1, NURBSPatch &p2)
Definition: nurbs.cpp:1073
void GetBdrElementTopo(Array< Element * > &boundary) const
Definition: nurbs.cpp:2136
static void Get3DRotationMatrix(double n[], double angle, double r, DenseMatrix &T)
Definition: nurbs.cpp:985
void SetKnotsFromPatches()
Definition: nurbs.cpp:2570
void GetElementTopo(Array< Element * > &elements) const
Definition: nurbs.cpp:2041
int GetNKV() const
Definition: nurbs.hpp:128
const Array< int > & GetOrders() const
Read-only access to the orders of all knot vectors.
Definition: nurbs.hpp:314
void UniformRefinement()
Definition: nurbs.cpp:541
Array< int > f_spaceOffsets
Definition: nurbs.hpp:194
int GetNBP() const
Definition: nurbs.hpp:311
Array< bool > activeBdrElem
Definition: nurbs.hpp:176
int GetNC() const
Definition: nurbs.hpp:127
void LoadSolution(std::istream &input, GridFunction &sol) const
Definition: nurbs.cpp:2607
int Dimension() const
Definition: mesh.hpp:645
KnotVector(const KnotVector &kv)
Definition: nurbs.hpp:44
int GetOrder() const
If all orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.
Definition: nurbs.hpp:317
KnotVector * GetKV(int i)
Definition: nurbs.hpp:129
int GetGNV() const
Definition: nurbs.hpp:321
void Print(std::ostream &out) const
Definition: nurbs.cpp:458
void LoadFE(int i, const FiniteElement *FE) const
Definition: nurbs.cpp:2509
Array< int > v_meshOffsets
Definition: nurbs.hpp:186
int GetNP() const
Definition: nurbs.hpp:310
void SetBdrPatchVertexMap(int p, const KnotVector *kv[], int *okv)
Definition: nurbs.cpp:3329
void GetBdrPatchKnotVectors(int p, Array< KnotVector * > &kv)
Definition: nurbs.cpp:1851
Array< int > p_meshOffsets
Definition: nurbs.hpp:189
KnotVector & operator=(const KnotVector &kv)
Definition: nurbs.cpp:45
void SetPatchVertexMap(int p, const KnotVector *kv[])
Definition: nurbs.cpp:3269
Table * GetElementDofTable()
Definition: nurbs.hpp:340
Table * GetBdrElementDofTable()
Definition: nurbs.hpp:341
Array2D< int > el_to_IJK
Definition: nurbs.hpp:201
void SetOrderFromOrders()
Definition: nurbs.cpp:1890
void CalcShape(Vector &shape, int i, double xi) const
Definition: nurbs.cpp:132
void CountBdrElements()
Definition: nurbs.cpp:2021
int GetNE() const
Definition: nurbs.hpp:324
void GetPatchNets(const Vector &Nodes, int vdim)
Definition: nurbs.cpp:2731
int MakeUniformDegree(int degree=-1)
Definition: nurbs.cpp:1050
bool isElement(int i) const
Definition: nurbs.hpp:57
Array< int > p_spaceOffsets
Definition: nurbs.hpp:195
void SetSolutionVector(Vector &Nodes, int vdim)
Definition: nurbs.cpp:2800
void Get3DElementTopo(Array< Element * > &elements) const
Definition: nurbs.cpp:2091
int findKnotSpan(double u) const
Definition: nurbs.cpp:214
void SetBdrPatchDofMap(int p, const KnotVector *kv[], int *okv)
Definition: nurbs.cpp:3358
int GetNCP() const
Definition: nurbs.hpp:50
double & operator()(int i, int j)
Definition: nurbs.hpp:451
void GenerateElementDofTable()
Definition: nurbs.cpp:2219
Array< NURBSPatch * > patches
Definition: nurbs.hpp:204
const double & operator[](int i) const
Definition: nurbs.hpp:81
void PrintCharacteristics(std::ostream &out) const
Definition: nurbs.cpp:1567
~KnotVector()
Destroys KnotVector.
Definition: nurbs.hpp:78
void Difference(const KnotVector &kv, Vector &diff) const
Definition: nurbs.cpp:243
Array< int > mOrders
Definition: nurbs.hpp:166
Vector data type.
Definition: vector.hpp:48
int operator[](const int i) const
Definition: nurbs.hpp:441
void Print(std::ostream &out) const
Definition: nurbs.cpp:125
double getKnotLocation(double xi, int ni) const
Definition: nurbs.hpp:59
int GetNBE() const
Definition: nurbs.hpp:326
int GetGNBE() const
Definition: nurbs.hpp:325
void CalcDShape(Vector &grad, int i, double xi) const
Definition: nurbs.cpp:158
void DegreeElevate(int rel_degree, int degree=16)
Definition: nurbs.cpp:2680
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
void CheckBdrPatches()
Definition: nurbs.cpp:1782
double & operator[](int i)
Definition: nurbs.hpp:80
void FlipDirection(int dir)
Definition: nurbs.cpp:948
void Get2DElementTopo(Array< Element * > &elements) const
Definition: nurbs.cpp:2055
GroupTopology gtopo
Definition: nurbs.hpp:386
int GetOrder() const
Definition: nurbs.hpp:51
Array< int > bel_to_patch
Definition: nurbs.hpp:200
void SetPatchDofMap(int p, const KnotVector *kv[])
Definition: nurbs.cpp:3299
void ConvertToPatches(const Vector &Nodes)
Definition: nurbs.cpp:2551
void KnotInsert(Array< KnotVector * > &kv)
Definition: nurbs.cpp:2704
void Print(std::ostream &out) const
Definition: nurbs.cpp:1532
void Get2DBdrElementTopo(Array< Element * > &boundary) const
Definition: nurbs.cpp:2150
int SetLoopDirection(int dir)
Definition: nurbs.cpp:482